Actual source code: itfunc.c

  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:   PCGetFailedReasonRank(pc,&pcreason);
216:   /* TODO: this code was wrong and is still wrong, there is no way to propogate the failure to all processes; their is no code to handle a ksp->reason on only some ranks */
217:   if (pcreason) {
218:     ksp->reason = KSP_DIVERGED_PC_FAILED;
219:   }
220:   return(0);
221: }

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

226:    Collective on ksp

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

232:    Level: intermediate

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

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

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

251:    Collective on ksp

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

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

259:    Level: intermediate

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

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

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

280:    Collective on ksp

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

286:    Level: intermediate

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

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

302:    Collective on ksp

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

307:    Level: developer

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


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

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

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

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

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

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

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

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

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

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

425: /*@C
426:    KSPConvergedReasonView - Displays the reason a KSP solve converged or diverged to a viewer

428:    Collective on ksp

430:    Parameter:
431: +  ksp - iterative context obtained from KSPCreate()
432: -  viewer - the viewer to display the reason

434:    Options Database Keys:
435: +  -ksp_converged_reason - print reason for converged or diverged, also prints number of iterations
436: -  -ksp_converged_reason ::failed - only print reason and number of iterations when diverged

438:    Notes:
439:      To change the format of the output call PetscViewerPushFormat(viewer,format) before this call. Use PETSC_VIEWER_DEFAULT for the default,
440:      use PETSC_VIEWER_FAILED to only display a reason if it fails.

442:    Level: beginner

444: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPConvergedDefault(),
445:           KSPSolveTranspose(), KSPGetIterationNumber(), KSP, KSPGetConvergedReason(), PetscViewerPushFormat(), PetscViewerPopFormat()
446: @*/
447: PetscErrorCode KSPConvergedReasonView(KSP ksp, PetscViewer viewer)
448: {
449:   PetscErrorCode    ierr;
450:   PetscBool         isAscii;
451:   PetscViewerFormat format;

454:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)ksp));
455:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
456:   if (isAscii) {
457:     PetscViewerGetFormat(viewer, &format);
458:     PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
459:     if (ksp->reason > 0 && format != PETSC_VIEWER_FAILED) {
460:       if (((PetscObject) ksp)->prefix) {
461:         PetscViewerASCIIPrintf(viewer,"Linear %s solve converged due to %s iterations %D\n",((PetscObject) ksp)->prefix,KSPConvergedReasons[ksp->reason],ksp->its);
462:       } else {
463:         PetscViewerASCIIPrintf(viewer,"Linear solve converged due to %s iterations %D\n",KSPConvergedReasons[ksp->reason],ksp->its);
464:       }
465:     } else if (ksp->reason <= 0) {
466:       if (((PetscObject) ksp)->prefix) {
467:         PetscViewerASCIIPrintf(viewer,"Linear %s solve did not converge due to %s iterations %D\n",((PetscObject) ksp)->prefix,KSPConvergedReasons[ksp->reason],ksp->its);
468:       } else {
469:         PetscViewerASCIIPrintf(viewer,"Linear solve did not converge due to %s iterations %D\n",KSPConvergedReasons[ksp->reason],ksp->its);
470:       }
471:       if (ksp->reason == KSP_DIVERGED_PC_FAILED) {
472:         PCFailedReason reason;
473:         PCGetFailedReason(ksp->pc,&reason);
474:         PetscViewerASCIIPrintf(viewer,"               PC failed due to %s \n",PCFailedReasons[reason]);
475:       }
476:     }
477:     PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
478:   }
479:   return(0);
480: }

482: /*@C
483:    KSPConvergedReasonViewSet - Sets an ADDITIONAL function that is to be used at the
484:     end of the linear solver to display the convergence reason of the linear solver.

486:    Logically Collective on KSP

488:    Input Parameters:
489: +  ksp - the KSP context
490: .  f - the ksp converged reason view function
491: .  vctx - [optional] user-defined context for private data for the
492:           ksp converged reason view routine (use NULL if no context is desired)
493: -  reasonviewdestroy - [optional] routine that frees reasonview context
494:           (may be NULL)

496:    Options Database Keys:
497: +    -ksp_converged_reason        - sets a default KSPConvergedReasonView()
498: -    -ksp_converged_reason_view_cancel - cancels all converged reason viewers that have
499:                             been hardwired into a code by
500:                             calls to KSPConvergedReasonViewSet(), but
501:                             does not cancel those set via
502:                             the options database.

504:    Notes:
505:    Several different converged reason view routines may be set by calling
506:    KSPConvergedReasonViewSet() multiple times; all will be called in the
507:    order in which they were set.

509:    Level: intermediate

511: .seealso: KSPConvergedReasonView(), KSPConvergedReasonViewCancel()
512: @*/
513: PetscErrorCode  KSPConvergedReasonViewSet(KSP ksp,PetscErrorCode (*f)(KSP,void*),void *vctx,PetscErrorCode (*reasonviewdestroy)(void**))
514: {
515:   PetscInt       i;
517:   PetscBool      identical;

521:   for (i=0; i<ksp->numberreasonviews;i++) {
522:     PetscMonitorCompare((PetscErrorCode (*)(void))f,vctx,reasonviewdestroy,(PetscErrorCode (*)(void))ksp->reasonview[i],ksp->reasonviewcontext[i],ksp->reasonviewdestroy[i],&identical);
523:     if (identical) return(0);
524:   }
525:   if (ksp->numberreasonviews >= MAXKSPREASONVIEWS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many KSP reasonview set");
526:   ksp->reasonview[ksp->numberreasonviews]          = f;
527:   ksp->reasonviewdestroy[ksp->numberreasonviews]   = reasonviewdestroy;
528:   ksp->reasonviewcontext[ksp->numberreasonviews++] = (void*)vctx;
529:   return(0);
530: }

532: /*@
533:    KSPConvergedReasonViewCancel - Clears all the reasonview functions for a KSP object.

535:    Collective on KSP

537:    Input Parameter:
538: .  ksp - iterative context obtained from KSPCreate()

540:    Level: intermediate

542: .seealso: KSPCreate(), KSPDestroy(), KSPReset()
543: @*/
544: PetscErrorCode  KSPConvergedReasonViewCancel(KSP ksp)
545: {
547:   PetscInt       i;

551:   for (i=0; i<ksp->numberreasonviews; i++) {
552:     if (ksp->reasonviewdestroy[i]) {
553:       (*ksp->reasonviewdestroy[i])(&ksp->reasonviewcontext[i]);
554:     }
555:   }
556:   ksp->numberreasonviews = 0;
557:   return(0);
558: }

560: /*@
561:   KSPConvergedReasonViewFromOptions - Processes command line options to determine if/how a KSPReason is to be viewed.

563:   Collective on ksp

565:   Input Parameters:
566: . ksp   - the KSP object

568:   Level: intermediate

570: .seealso: KSPConvergedReasonView()
571: @*/
572: PetscErrorCode KSPConvergedReasonViewFromOptions(KSP ksp)
573: {
574:   PetscViewer       viewer;
575:   PetscBool         flg;
576:   PetscViewerFormat format;
577:   PetscErrorCode    ierr;
578:   PetscInt          i;


582:   /* Call all user-provided reason review routines */
583:   for (i=0; i<ksp->numberreasonviews; i++) {
584:     (*ksp->reasonview[i])(ksp,ksp->reasonviewcontext[i]);
585:   }

587:   /* Call the default PETSc routine */
588:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->options,((PetscObject)ksp)->prefix,"-ksp_converged_reason",&viewer,&format,&flg);
589:   if (flg) {
590:     PetscViewerPushFormat(viewer,format);
591:     KSPConvergedReasonView(ksp, viewer);
592:     PetscViewerPopFormat(viewer);
593:     PetscViewerDestroy(&viewer);
594:   }
595:   return(0);
596: }

598: /*@C
599:   KSPConvergedRateView - Displays the reason a KSP solve converged or diverged to a viewer

601:   Collective on ksp

603:   Input Parameters:
604: +  ksp    - iterative context obtained from KSPCreate()
605: -  viewer - the viewer to display the reason

607:   Options Database Keys:
608: . -ksp_converged_rate - print reason for convergence or divergence and the convergence rate (or 0.0 for divergence)

610:   Notes:
611:   To change the format of the output, call PetscViewerPushFormat(viewer,format) before this call.

613:   Suppose that the residual is reduced linearly, $r_k = c^k r_0$, which means $log r_k = log r_0 + k log c$. After linear regression,
614:   the slope is $\log c$. The coefficient of determination is given by $1 - \frac{\sum_i (y_i - f(x_i))^2}{\sum_i (y_i - \bar y)}$,
615:   see also https://en.wikipedia.org/wiki/Coefficient_of_determination

617:   Level: intermediate

619: .seealso: KSPConvergedReasonView(), KSPGetConvergedRate(), KSPSetTolerances(), KSPConvergedDefault()
620: @*/
621: PetscErrorCode KSPConvergedRateView(KSP ksp, PetscViewer viewer)
622: {
623:   PetscViewerFormat format;
624:   PetscBool         isAscii;
625:   PetscReal         rrate, rRsq, erate = 0.0, eRsq = 0.0;
626:   PetscInt          its;
627:   const char       *prefix, *reason = KSPConvergedReasons[ksp->reason];
628:   PetscErrorCode    ierr;

631:   KSPGetOptionsPrefix(ksp, &prefix);
632:   KSPGetIterationNumber(ksp, &its);
633:   KSPComputeConvergenceRate(ksp, &rrate, &rRsq, &erate, &eRsq);
634:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)ksp));
635:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
636:   if (isAscii) {
637:     PetscViewerGetFormat(viewer, &format);
638:     PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
639:     if (ksp->reason > 0) {
640:       if (prefix) {PetscViewerASCIIPrintf(viewer, "Linear %s solve converged due to %s iterations %D", prefix, reason, its);}
641:       else        {PetscViewerASCIIPrintf(viewer, "Linear solve converged due to %s iterations %D", reason, its);}
642:       PetscViewerASCIIUseTabs(viewer, PETSC_FALSE);
643:       if (rRsq >= 0.0) {PetscViewerASCIIPrintf(viewer, " res rate %g R^2 %g", rrate, rRsq);}
644:       if (eRsq >= 0.0) {PetscViewerASCIIPrintf(viewer, " error rate %g R^2 %g", erate, eRsq);}
645:       PetscViewerASCIIPrintf(viewer, "\n");
646:       PetscViewerASCIIUseTabs(viewer, PETSC_TRUE);
647:     } else if (ksp->reason <= 0) {
648:       if (prefix) {PetscViewerASCIIPrintf(viewer, "Linear %s solve did not converge due to %s iterations %D", prefix, reason, its);}
649:       else        {PetscViewerASCIIPrintf(viewer, "Linear solve did not converge due to %s iterations %D", reason, its);}
650:       PetscViewerASCIIUseTabs(viewer, PETSC_FALSE);
651:       if (rRsq >= 0.0) {PetscViewerASCIIPrintf(viewer, " res rate %g R^2 %g", rrate, rRsq);}
652:       if (eRsq >= 0.0) {PetscViewerASCIIPrintf(viewer, " error rate %g R^2 %g", erate, eRsq);}
653:       PetscViewerASCIIPrintf(viewer, "\n");
654:       PetscViewerASCIIUseTabs(viewer, PETSC_TRUE);
655:       if (ksp->reason == KSP_DIVERGED_PC_FAILED) {
656:         PCFailedReason reason;
657:         PCGetFailedReason(ksp->pc,&reason);
658:         PetscViewerASCIIPrintf(viewer,"               PC failed due to %s \n",PCFailedReasons[reason]);
659:       }
660:     }
661:     PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
662:   }
663:   return(0);
664: }

666: #include <petscdraw.h>

668: static PetscErrorCode KSPViewEigenvalues_Internal(KSP ksp, PetscBool isExplicit, PetscViewer viewer, PetscViewerFormat format)
669: {
670:   PetscReal     *r, *c;
671:   PetscInt       n, i, neig;
672:   PetscBool      isascii, isdraw;
673:   PetscMPIInt    rank;

677:   MPI_Comm_rank(PetscObjectComm((PetscObject) ksp), &rank);
678:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);
679:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERDRAW,  &isdraw);
680:   if (isExplicit) {
681:     VecGetSize(ksp->vec_sol,&n);
682:     PetscMalloc2(n, &r, n, &c);
683:     KSPComputeEigenvaluesExplicitly(ksp, n, r, c);
684:     neig = n;
685:   } else {
686:     PetscInt nits;

688:     KSPGetIterationNumber(ksp, &nits);
689:     n    = nits+2;
690:     if (!nits) {PetscViewerASCIIPrintf(viewer, "Zero iterations in solver, cannot approximate any eigenvalues\n");return(0);}
691:     PetscMalloc2(n, &r, n, &c);
692:     KSPComputeEigenvalues(ksp, n, r, c, &neig);
693:   }
694:   if (isascii) {
695:     PetscViewerASCIIPrintf(viewer, "%s computed eigenvalues\n", isExplicit ? "Explicitly" : "Iteratively");
696:     for (i = 0; i < neig; ++i) {
697:       if (c[i] >= 0.0) {PetscViewerASCIIPrintf(viewer, "%g + %gi\n", (double) r[i],  (double) c[i]);}
698:       else             {PetscViewerASCIIPrintf(viewer, "%g - %gi\n", (double) r[i], -(double) c[i]);}
699:     }
700:   } else if (isdraw && !rank) {
701:     PetscDraw   draw;
702:     PetscDrawSP drawsp;

704:     if (format == PETSC_VIEWER_DRAW_CONTOUR) {
705:       KSPPlotEigenContours_Private(ksp,neig,r,c);
706:     } else {
707:       if (!ksp->eigviewer) {PetscViewerDrawOpen(PETSC_COMM_SELF,NULL,isExplicit ? "Explicitly Computed Eigenvalues" : "Iteratively Computed Eigenvalues",PETSC_DECIDE,PETSC_DECIDE,400,400,&ksp->eigviewer);}
708:       PetscViewerDrawGetDraw(ksp->eigviewer,0,&draw);
709:       PetscDrawSPCreate(draw,1,&drawsp);
710:       PetscDrawSPReset(drawsp);
711:       for (i = 0; i < neig; ++i) {PetscDrawSPAddPoint(drawsp,r+i,c+i);}
712:       PetscDrawSPDraw(drawsp,PETSC_TRUE);
713:       PetscDrawSPSave(drawsp);
714:       PetscDrawSPDestroy(&drawsp);
715:     }
716:   }
717:   PetscFree2(r, c);
718:   return(0);
719: }

721: static PetscErrorCode KSPViewSingularvalues_Internal(KSP ksp, PetscViewer viewer, PetscViewerFormat format)
722: {
723:   PetscReal      smax, smin;
724:   PetscInt       nits;
725:   PetscBool      isascii;

729:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);
730:   KSPGetIterationNumber(ksp, &nits);
731:   if (!nits) {PetscViewerASCIIPrintf(viewer, "Zero iterations in solver, cannot approximate any singular values\n");return(0);}
732:   KSPComputeExtremeSingularValues(ksp, &smax, &smin);
733:   if (isascii) {PetscViewerASCIIPrintf(viewer, "Iteratively computed extreme singular values: max %g min %g max/min %g\n",(double)smax,(double)smin,(double)(smax/smin));}
734:   return(0);
735: }

737: static PetscErrorCode KSPViewFinalResidual_Internal(KSP ksp, PetscViewer viewer, PetscViewerFormat format)
738: {
739:   PetscBool      isascii;

743:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);
744:   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");
745:   if (isascii) {
746:     Mat       A;
747:     Vec       t;
748:     PetscReal norm;

750:     PCGetOperators(ksp->pc, &A, NULL);
751:     VecDuplicate(ksp->vec_rhs, &t);
752:     KSP_MatMult(ksp, A, ksp->vec_sol, t);
753:     VecAYPX(t, -1.0, ksp->vec_rhs);
754:     VecNorm(t, NORM_2, &norm);
755:     VecDestroy(&t);
756:     PetscViewerASCIIPrintf(viewer, "KSP final norm of residual %g\n", (double) norm);
757:   }
758:   return(0);
759: }

761: static PetscErrorCode KSPMonitorPauseFinal_Internal(KSP ksp)
762: {
763:   PetscInt       i;

767:   if (!ksp->pauseFinal) return(0);
768:   for (i = 0; i < ksp->numbermonitors; ++i) {
769:     PetscViewerAndFormat *vf = (PetscViewerAndFormat *) ksp->monitorcontext[i];
770:     PetscDraw             draw;
771:     PetscReal             lpause;

773:     if (!vf) continue;
774:     if (vf->lg) {
776:       if (((PetscObject) vf->lg)->classid != PETSC_DRAWLG_CLASSID) continue;
777:       PetscDrawLGGetDraw(vf->lg, &draw);
778:       PetscDrawGetPause(draw, &lpause);
779:       PetscDrawSetPause(draw, -1.0);
780:       PetscDrawPause(draw);
781:       PetscDrawSetPause(draw, lpause);
782:     } else {
783:       PetscBool isdraw;

786:       if (((PetscObject) vf->viewer)->classid != PETSC_VIEWER_CLASSID) continue;
787:       PetscObjectTypeCompare((PetscObject) vf->viewer, PETSCVIEWERDRAW, &isdraw);
788:       if (!isdraw) continue;
789:       PetscViewerDrawGetDraw(vf->viewer, 0, &draw);
790:       PetscDrawGetPause(draw, &lpause);
791:       PetscDrawSetPause(draw, -1.0);
792:       PetscDrawPause(draw);
793:       PetscDrawSetPause(draw, lpause);
794:     }
795:   }
796:   return(0);
797: }

799: static PetscErrorCode KSPSolve_Private(KSP ksp,Vec b,Vec x)
800: {
802:   PetscBool      flg = PETSC_FALSE,inXisinB=PETSC_FALSE,guess_zero;
803:   Mat            mat,pmat;
804:   MPI_Comm       comm;
805:   MatNullSpace   nullsp;
806:   Vec            btmp,vec_rhs=NULL;

809:   comm = PetscObjectComm((PetscObject)ksp);
810:   if (x && x == b) {
811:     if (!ksp->guess_zero) SETERRQ(comm,PETSC_ERR_ARG_INCOMP,"Cannot use x == b with nonzero initial guess");
812:     VecDuplicate(b,&x);
813:     inXisinB = PETSC_TRUE;
814:   }
815:   if (b) {
816:     PetscObjectReference((PetscObject)b);
817:     VecDestroy(&ksp->vec_rhs);
818:     ksp->vec_rhs = b;
819:   }
820:   if (x) {
821:     PetscObjectReference((PetscObject)x);
822:     VecDestroy(&ksp->vec_sol);
823:     ksp->vec_sol = x;
824:   }

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

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

830:   /* reset the residual history list if requested */
831:   if (ksp->res_hist_reset) ksp->res_hist_len = 0;
832:   if (ksp->err_hist_reset) ksp->err_hist_len = 0;

834:   if (ksp->guess) {
835:     PetscObjectState ostate,state;

837:     KSPGuessSetUp(ksp->guess);
838:     PetscObjectStateGet((PetscObject)ksp->vec_sol,&ostate);
839:     KSPGuessFormGuess(ksp->guess,ksp->vec_rhs,ksp->vec_sol);
840:     PetscObjectStateGet((PetscObject)ksp->vec_sol,&state);
841:     if (state != ostate) {
842:       ksp->guess_zero = PETSC_FALSE;
843:     } else {
844:       PetscInfo(ksp,"Using zero initial guess since the KSPGuess object did not change the vector\n");
845:       ksp->guess_zero = PETSC_TRUE;
846:     }
847:   }

849:   /* KSPSetUp() scales the matrix if needed */
850:   KSPSetUp(ksp);
851:   KSPSetUpOnBlocks(ksp);

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

855:   PetscLogEventBegin(KSP_Solve,ksp,ksp->vec_rhs,ksp->vec_sol,0);
856:   PCGetOperators(ksp->pc,&mat,&pmat);
857:   /* diagonal scale RHS if called for */
858:   if (ksp->dscale) {
859:     VecPointwiseMult(ksp->vec_rhs,ksp->vec_rhs,ksp->diagonal);
860:     /* second time in, but matrix was scaled back to original */
861:     if (ksp->dscalefix && ksp->dscalefix2) {
862:       Mat mat,pmat;

864:       PCGetOperators(ksp->pc,&mat,&pmat);
865:       MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
866:       if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
867:     }

869:     /* scale initial guess */
870:     if (!ksp->guess_zero) {
871:       if (!ksp->truediagonal) {
872:         VecDuplicate(ksp->diagonal,&ksp->truediagonal);
873:         VecCopy(ksp->diagonal,ksp->truediagonal);
874:         VecReciprocal(ksp->truediagonal);
875:       }
876:       VecPointwiseMult(ksp->vec_sol,ksp->vec_sol,ksp->truediagonal);
877:     }
878:   }
879:   PCPreSolve(ksp->pc,ksp);

881:   if (ksp->guess_zero) { VecSet(ksp->vec_sol,0.0);}
882:   if (ksp->guess_knoll) { /* The Knoll trick is independent on the KSPGuess specified */
883:     PCApply(ksp->pc,ksp->vec_rhs,ksp->vec_sol);
884:     KSP_RemoveNullSpace(ksp,ksp->vec_sol);
885:     ksp->guess_zero = PETSC_FALSE;
886:   }

888:   /* can we mark the initial guess as zero for this solve? */
889:   guess_zero = ksp->guess_zero;
890:   if (!ksp->guess_zero) {
891:     PetscReal norm;

893:     VecNormAvailable(ksp->vec_sol,NORM_2,&flg,&norm);
894:     if (flg && !norm) ksp->guess_zero = PETSC_TRUE;
895:   }
896:   if (ksp->transpose_solve) {
897:     MatGetNullSpace(pmat,&nullsp);
898:   } else {
899:     MatGetTransposeNullSpace(pmat,&nullsp);
900:   }
901:   if (nullsp) {
902:     VecDuplicate(ksp->vec_rhs,&btmp);
903:     VecCopy(ksp->vec_rhs,btmp);
904:     MatNullSpaceRemove(nullsp,btmp);
905:     vec_rhs      = ksp->vec_rhs;
906:     ksp->vec_rhs = btmp;
907:   }
908:   VecLockReadPush(ksp->vec_rhs);
909:   if (ksp->reason == KSP_DIVERGED_PC_FAILED) {
910:     VecSetInf(ksp->vec_sol);
911:   }
912:   (*ksp->ops->solve)(ksp);
913:   KSPMonitorPauseFinal_Internal(ksp);

915:   VecLockReadPop(ksp->vec_rhs);
916:   if (nullsp) {
917:     ksp->vec_rhs = vec_rhs;
918:     VecDestroy(&btmp);
919:   }

921:   ksp->guess_zero = guess_zero;

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

926:   KSPConvergedReasonViewFromOptions(ksp);

928:   if (ksp->viewRate) {
929:     PetscViewerPushFormat(ksp->viewerRate,ksp->formatRate);
930:     KSPConvergedRateView(ksp, ksp->viewerRate);
931:     PetscViewerPopFormat(ksp->viewerRate);
932:   }
933:   PCPostSolve(ksp->pc,ksp);

935:   /* diagonal scale solution if called for */
936:   if (ksp->dscale) {
937:     VecPointwiseMult(ksp->vec_sol,ksp->vec_sol,ksp->diagonal);
938:     /* unscale right hand side and matrix */
939:     if (ksp->dscalefix) {
940:       Mat mat,pmat;

942:       VecReciprocal(ksp->diagonal);
943:       VecPointwiseMult(ksp->vec_rhs,ksp->vec_rhs,ksp->diagonal);
944:       PCGetOperators(ksp->pc,&mat,&pmat);
945:       MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
946:       if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
947:       VecReciprocal(ksp->diagonal);
948:       ksp->dscalefix2 = PETSC_TRUE;
949:     }
950:   }
951:   PetscLogEventEnd(KSP_Solve,ksp,ksp->vec_rhs,ksp->vec_sol,0);
952:   if (ksp->guess) {
953:     KSPGuessUpdate(ksp->guess,ksp->vec_rhs,ksp->vec_sol);
954:   }
955:   if (ksp->postsolve) {
956:     (*ksp->postsolve)(ksp,ksp->vec_rhs,ksp->vec_sol,ksp->postctx);
957:   }

959:   PCGetOperators(ksp->pc,&mat,&pmat);
960:   if (ksp->viewEV)       {KSPViewEigenvalues_Internal(ksp, PETSC_FALSE, ksp->viewerEV,    ksp->formatEV);}
961:   if (ksp->viewEVExp)    {KSPViewEigenvalues_Internal(ksp, PETSC_TRUE,  ksp->viewerEVExp, ksp->formatEVExp);}
962:   if (ksp->viewSV)       {KSPViewSingularvalues_Internal(ksp, ksp->viewerSV, ksp->formatSV);}
963:   if (ksp->viewFinalRes) {KSPViewFinalResidual_Internal(ksp, ksp->viewerFinalRes, ksp->formatFinalRes);}
964:   if (ksp->viewMat)      {ObjectView((PetscObject) mat,           ksp->viewerMat,    ksp->formatMat);}
965:   if (ksp->viewPMat)     {ObjectView((PetscObject) pmat,          ksp->viewerPMat,   ksp->formatPMat);}
966:   if (ksp->viewRhs)      {ObjectView((PetscObject) ksp->vec_rhs,  ksp->viewerRhs,    ksp->formatRhs);}
967:   if (ksp->viewSol)      {ObjectView((PetscObject) ksp->vec_sol,  ksp->viewerSol,    ksp->formatSol);}
968:   if (ksp->view)         {ObjectView((PetscObject) ksp,           ksp->viewer,       ksp->format);}
969:   if (ksp->viewDScale)   {ObjectView((PetscObject) ksp->diagonal, ksp->viewerDScale, ksp->formatDScale);}
970:   if (ksp->viewMatExp)   {
971:     Mat A, B;

973:     PCGetOperators(ksp->pc, &A, NULL);
974:     if (ksp->transpose_solve) {
975:       Mat AT;

977:       MatCreateTranspose(A, &AT);
978:       MatComputeOperator(AT, MATAIJ, &B);
979:       MatDestroy(&AT);
980:     } else {
981:       MatComputeOperator(A, MATAIJ, &B);
982:     }
983:     ObjectView((PetscObject) B, ksp->viewerMatExp, ksp->formatMatExp);
984:     MatDestroy(&B);
985:   }
986:   if (ksp->viewPOpExp)   {
987:     Mat B;

989:     KSPComputeOperator(ksp, MATAIJ, &B);
990:     ObjectView((PetscObject) B, ksp->viewerPOpExp, ksp->formatPOpExp);
991:     MatDestroy(&B);
992:   }

994:   if (inXisinB) {
995:     VecCopy(x,b);
996:     VecDestroy(&x);
997:   }
998:   PetscObjectSAWsBlock((PetscObject)ksp);
999:   if (ksp->errorifnotconverged && ksp->reason < 0 && ksp->reason != KSP_DIVERGED_ITS) {
1000:     if (ksp->reason == KSP_DIVERGED_PC_FAILED) {
1001:       PCFailedReason reason;
1002:       PCGetFailedReason(ksp->pc,&reason);
1003:       SETERRQ2(comm,PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged, reason %s PC failed due to %s",KSPConvergedReasons[ksp->reason],PCFailedReasons[reason]);
1004:     } else SETERRQ1(comm,PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged, reason %s",KSPConvergedReasons[ksp->reason]);
1005:   }
1006:   return(0);
1007: }

1009: /*@
1010:    KSPSolve - Solves linear system.

1012:    Collective on ksp

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

1019:    Options Database Keys:
1020: +  -ksp_view_eigenvalues - compute preconditioned operators eigenvalues
1021: .  -ksp_view_eigenvalues_explicit - compute the eigenvalues by forming the dense operator and using LAPACK
1022: .  -ksp_view_mat binary - save matrix to the default binary viewer
1023: .  -ksp_view_pmat binary - save matrix used to build preconditioner to the default binary viewer
1024: .  -ksp_view_rhs binary - save right hand side vector to the default binary viewer
1025: .  -ksp_view_solution binary - save computed solution vector to the default binary viewer
1026:            (can be read later with src/ksp/tutorials/ex10.c for testing solvers)
1027: .  -ksp_view_mat_explicit - for matrix-free operators, computes the matrix entries and views them
1028: .  -ksp_view_preconditioned_operator_explicit - computes the product of the preconditioner and matrix as an explicit matrix and views it
1029: .  -ksp_converged_reason - print reason for converged or diverged, also prints number of iterations
1030: .  -ksp_view_final_residual - print 2-norm of true linear system residual at the end of the solution process
1031: -  -ksp_view - print the ksp data structure at the end of the system solution

1033:    Notes:

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

1037:    The operator is specified with KSPSetOperators().

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

1042:    If you provide a matrix that has a MatSetNullSpace() and MatSetTransposeNullSpace() this will use that information to solve singular systems
1043:    in the least squares sense with a norm minimizing solution.
1044: $
1045: $                   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()
1046: $
1047: $    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
1048: $    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
1049: $    direction thus the solution which is a linear combination of the search directions has no component in the nullspace(A).
1050: $
1051: $    We recommend always using GMRES for such singular systems.
1052: $    If nullspace(A) = nullspace(A') (note symmetric matrices always satisfy this property) then both left and right preconditioning will work
1053: $    If nullspace(A) != nullspace(A') then left preconditioning will work but right preconditioning may not work (or it may).

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


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

1065:    Understanding Convergence:
1066:    The routines KSPMonitorSet(), KSPComputeEigenvalues(), and
1067:    KSPComputeEigenvaluesExplicitly() provide information on additional
1068:    options to monitor convergence and print eigenvalue information.

1070:    Level: beginner

1072: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPConvergedDefault(),
1073:           KSPSolveTranspose(), KSPGetIterationNumber(), MatNullSpaceCreate(), MatSetNullSpace(), MatSetTransposeNullSpace(), KSP,
1074:           KSPConvergedReasonView()
1075: @*/
1076: PetscErrorCode KSPSolve(KSP ksp,Vec b,Vec x)
1077: {

1084:   ksp->transpose_solve = PETSC_FALSE;
1085:   KSPSolve_Private(ksp,b,x);
1086:   return(0);
1087: }

1089: /*@
1090:    KSPSolveTranspose - Solves the transpose of a linear system.

1092:    Collective on ksp

1094:    Input Parameters:
1095: +  ksp - iterative context obtained from KSPCreate()
1096: .  b - right hand side vector
1097: -  x - solution vector

1099:    Notes:
1100:     For complex numbers this solve the non-Hermitian transpose system.

1102:    Developer Notes:
1103:     We need to implement a KSPSolveHermitianTranspose()

1105:    Level: developer

1107: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPConvergedDefault(),
1108:           KSPSolve(), KSP
1109: @*/
1110: PetscErrorCode KSPSolveTranspose(KSP ksp,Vec b,Vec x)
1111: {

1118:   if (ksp->transpose.use_explicittranspose) {
1119:     Mat J,Jpre;
1120:     KSPGetOperators(ksp,&J,&Jpre);
1121:     if (!ksp->transpose.reuse_transpose) {
1122:       MatTranspose(J,MAT_INITIAL_MATRIX,&ksp->transpose.AT);
1123:       if (J != Jpre) {
1124:         MatTranspose(Jpre,MAT_INITIAL_MATRIX,&ksp->transpose.BT);
1125:       }
1126:       ksp->transpose.reuse_transpose = PETSC_TRUE;
1127:     } else {
1128:       MatTranspose(J,MAT_REUSE_MATRIX,&ksp->transpose.AT);
1129:       if (J != Jpre) {
1130:         MatTranspose(Jpre,MAT_REUSE_MATRIX,&ksp->transpose.BT);
1131:       }
1132:     }
1133:     if (J == Jpre && ksp->transpose.BT != ksp->transpose.AT) {
1134:       PetscObjectReference((PetscObject)ksp->transpose.AT);
1135:       ksp->transpose.BT = ksp->transpose.AT;
1136:     }
1137:     KSPSetOperators(ksp,ksp->transpose.AT,ksp->transpose.BT);
1138:   } else {
1139:     ksp->transpose_solve = PETSC_TRUE;
1140:   }
1141:   KSPSolve_Private(ksp,b,x);
1142:   return(0);
1143: }

1145: static PetscErrorCode KSPViewFinalMatResidual_Internal(KSP ksp, Mat B, Mat X, PetscViewer viewer, PetscViewerFormat format, PetscInt shift)
1146: {
1147:   Mat            A, R;
1148:   PetscReal      *norms;
1149:   PetscInt       i, N;
1150:   PetscBool      flg;

1154:   PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &flg);
1155:   if (flg) {
1156:     PCGetOperators(ksp->pc, &A, NULL);
1157:     MatMatMult(A, X, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &R);
1158:     MatAYPX(R, -1.0, B, SAME_NONZERO_PATTERN);
1159:     MatGetSize(R, NULL, &N);
1160:     PetscMalloc1(N, &norms);
1161:     MatGetColumnNorms(R, NORM_2, norms);
1162:     MatDestroy(&R);
1163:     for (i = 0; i < N; ++i) {
1164:       PetscViewerASCIIPrintf(viewer, "%s #%D %g\n", i == 0 ? "KSP final norm of residual" : "                          ", shift + i, (double)norms[i]);
1165:     }
1166:     PetscFree(norms);
1167:   }
1168:   return(0);
1169: }

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

1174:    Input Parameters:
1175: +     ksp - iterative context
1176: -     B - block of right-hand sides

1178:    Output Parameter:
1179: .     X - block of solutions

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

1184:    Level: intermediate

1186: .seealso:  KSPSolve(), MatMatSolve(), MATDENSE, KSPHPDDM, PCBJACOBI, PCASM
1187: @*/
1188: PetscErrorCode KSPMatSolve(KSP ksp, Mat B, Mat X)
1189: {
1190:   Mat            A, vB, vX;
1191:   Vec            cb, cx;
1192:   PetscInt       n1, N1, n2, N2, Bbn = PETSC_DECIDE;
1193:   PetscBool      match;

1202:   if (!B->assembled) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
1203:   MatCheckPreallocated(X, 3);
1204:   if (!X->assembled) {
1205:     MatSetOption(X, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE);
1206:     MatAssemblyBegin(X, MAT_FINAL_ASSEMBLY);
1207:     MatAssemblyEnd(X, MAT_FINAL_ASSEMBLY);
1208:   }
1209:   if (B == X) SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_IDN, "B and X must be different matrices");
1210:   KSPGetOperators(ksp, &A, NULL);
1211:   MatGetLocalSize(B, NULL, &n2);
1212:   MatGetLocalSize(X, NULL, &n1);
1213:   MatGetSize(B, NULL, &N2);
1214:   MatGetSize(X, NULL, &N1);
1215:   if (n1 != n2 || N1 != N2) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible number of columns between block of right-hand sides (n,N) = (%D,%D) and block of solutions (n,N) = (%D,%D)", n2, N2, n1, N1);
1216:   PetscObjectBaseTypeCompareAny((PetscObject)B, &match, MATSEQDENSE, MATMPIDENSE, "");
1217:   if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of right-hand sides not stored in a dense Mat");
1218:   PetscObjectBaseTypeCompareAny((PetscObject)X, &match, MATSEQDENSE, MATMPIDENSE, "");
1219:   if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of solutions not stored in a dense Mat");
1220:   KSPSetUp(ksp);
1221:   KSPSetUpOnBlocks(ksp);
1222:   if (ksp->ops->matsolve) {
1223:     if (ksp->guess_zero) {
1224:       MatZeroEntries(X);
1225:     }
1226:     PetscLogEventBegin(KSP_MatSolve, ksp, B, X, 0);
1227:     KSPGetMatSolveBatchSize(ksp, &Bbn);
1228:     /* by default, do a single solve with all columns */
1229:     if (Bbn == PETSC_DECIDE) Bbn = N2;
1230:     else if (Bbn < 1) SETERRQ1(PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "KSPMatSolve() block size %D must be positive", Bbn);
1231:     PetscInfo2(ksp, "KSP type %s solving using blocks of width at most %D\n", ((PetscObject)ksp)->type_name, Bbn);
1232:     /* if -ksp_matsolve_batch_size is greater than the actual number of columns, do a single solve with all columns */
1233:     if (Bbn >= N2) {
1234:       (*ksp->ops->matsolve)(ksp, B, X);
1235:       if (ksp->viewFinalRes) {
1236:         KSPViewFinalMatResidual_Internal(ksp, B, X, ksp->viewerFinalRes, ksp->formatFinalRes, 0);
1237:       }

1239:       KSPConvergedReasonViewFromOptions(ksp);

1241:       if (ksp->viewRate) {
1242:         PetscViewerPushFormat(ksp->viewerRate,PETSC_VIEWER_DEFAULT);
1243:         KSPConvergedRateView(ksp, ksp->viewerRate);
1244:         PetscViewerPopFormat(ksp->viewerRate);
1245:       }
1246:     } else {
1247:       for (n2 = 0; n2 < N2; n2 += Bbn) {
1248:         MatDenseGetSubMatrix(B, n2, PetscMin(n2+Bbn, N2), &vB);
1249:         MatDenseGetSubMatrix(X, n2, PetscMin(n2+Bbn, N2), &vX);
1250:         (*ksp->ops->matsolve)(ksp, vB, vX);
1251:         if (ksp->viewFinalRes) {
1252:           KSPViewFinalMatResidual_Internal(ksp, vB, vX, ksp->viewerFinalRes, ksp->formatFinalRes, n2);
1253:         }

1255:         KSPConvergedReasonViewFromOptions(ksp);

1257:         if (ksp->viewRate) {
1258:           PetscViewerPushFormat(ksp->viewerRate,PETSC_VIEWER_DEFAULT);
1259:           KSPConvergedRateView(ksp, ksp->viewerRate);
1260:           PetscViewerPopFormat(ksp->viewerRate);
1261:         }
1262:         MatDenseRestoreSubMatrix(B, &vB);
1263:         MatDenseRestoreSubMatrix(X, &vX);
1264:       }
1265:     }
1266:     if (ksp->view) {
1267:       KSPView(ksp, ksp->viewer);
1268:     }
1269:     PetscLogEventEnd(KSP_MatSolve, ksp, B, X, 0);
1270:   } else {
1271:     PetscInfo1(ksp, "KSP type %s solving column by column\n", ((PetscObject)ksp)->type_name);
1272:     for (n2 = 0; n2 < N2; ++n2) {
1273:       MatDenseGetColumnVecRead(B, n2, &cb);
1274:       MatDenseGetColumnVecWrite(X, n2, &cx);
1275:       KSPSolve(ksp, cb, cx);
1276:       MatDenseRestoreColumnVecWrite(X, n2, &cx);
1277:       MatDenseRestoreColumnVecRead(B, n2, &cb);
1278:     }
1279:   }
1280:   return(0);
1281: }

1283: /*@
1284:      KSPSetMatSolveBatchSize - Sets the maximum number of columns treated simultaneously in KSPMatSolve().

1286:     Logically collective

1288:    Input Parameters:
1289: +     ksp - iterative context
1290: -     bs - block size

1292:    Level: advanced

1294: .seealso:  KSPMatSolve(), KSPGetMatSolveBatchSize(), -mat_mumps_icntl_27, -matmatmult_Bbn
1295: @*/
1296: PetscErrorCode KSPSetMatSolveBatchSize(KSP ksp, PetscInt bs)
1297: {
1301:   ksp->nmax = bs;
1302:   return(0);
1303: }

1305: /*@
1306:      KSPGetMatSolveBatchSize - Gets the maximum number of columns treated simultaneously in KSPMatSolve().

1308:    Input Parameter:
1309: .     ksp - iterative context

1311:    Output Parameter:
1312: .     bs - block size

1314:    Level: advanced

1316: .seealso:  KSPMatSolve(), KSPSetMatSolveBatchSize(), -mat_mumps_icntl_27, -matmatmult_Bbn
1317: @*/
1318: PetscErrorCode KSPGetMatSolveBatchSize(KSP ksp, PetscInt *bs)
1319: {
1323:   *bs = ksp->nmax;
1324:   return(0);
1325: }

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

1330:    Collective on ksp

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

1335:    Level: beginner

1337: .seealso: KSPCreate(), KSPSetUp(), KSPSolve(), KSPSetFromOptions(), KSP
1338: @*/
1339: PetscErrorCode  KSPResetViewers(KSP ksp)
1340: {

1345:   if (!ksp) return(0);
1346:   PetscViewerDestroy(&ksp->viewer);
1347:   PetscViewerDestroy(&ksp->viewerPre);
1348:   PetscViewerDestroy(&ksp->viewerRate);
1349:   PetscViewerDestroy(&ksp->viewerMat);
1350:   PetscViewerDestroy(&ksp->viewerPMat);
1351:   PetscViewerDestroy(&ksp->viewerRhs);
1352:   PetscViewerDestroy(&ksp->viewerSol);
1353:   PetscViewerDestroy(&ksp->viewerMatExp);
1354:   PetscViewerDestroy(&ksp->viewerEV);
1355:   PetscViewerDestroy(&ksp->viewerSV);
1356:   PetscViewerDestroy(&ksp->viewerEVExp);
1357:   PetscViewerDestroy(&ksp->viewerFinalRes);
1358:   PetscViewerDestroy(&ksp->viewerPOpExp);
1359:   PetscViewerDestroy(&ksp->viewerDScale);
1360:   ksp->view         = PETSC_FALSE;
1361:   ksp->viewPre      = PETSC_FALSE;
1362:   ksp->viewMat      = PETSC_FALSE;
1363:   ksp->viewPMat     = PETSC_FALSE;
1364:   ksp->viewRhs      = PETSC_FALSE;
1365:   ksp->viewSol      = PETSC_FALSE;
1366:   ksp->viewMatExp   = PETSC_FALSE;
1367:   ksp->viewEV       = PETSC_FALSE;
1368:   ksp->viewSV       = PETSC_FALSE;
1369:   ksp->viewEVExp    = PETSC_FALSE;
1370:   ksp->viewFinalRes = PETSC_FALSE;
1371:   ksp->viewPOpExp   = PETSC_FALSE;
1372:   ksp->viewDScale   = PETSC_FALSE;
1373:   return(0);
1374: }

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

1379:    Collective on ksp

1381:    Input Parameter:
1382: .  ksp - iterative context obtained from KSPCreate()

1384:    Level: beginner

1386: .seealso: KSPCreate(), KSPSetUp(), KSPSolve(), KSP
1387: @*/
1388: PetscErrorCode  KSPReset(KSP ksp)
1389: {

1394:   if (!ksp) return(0);
1395:   if (ksp->ops->reset) {
1396:     (*ksp->ops->reset)(ksp);
1397:   }
1398:   if (ksp->pc) {PCReset(ksp->pc);}
1399:   if (ksp->guess) {
1400:     KSPGuess guess = ksp->guess;
1401:     if (guess->ops->reset) { (*guess->ops->reset)(guess); }
1402:   }
1403:   VecDestroyVecs(ksp->nwork,&ksp->work);
1404:   VecDestroy(&ksp->vec_rhs);
1405:   VecDestroy(&ksp->vec_sol);
1406:   VecDestroy(&ksp->diagonal);
1407:   VecDestroy(&ksp->truediagonal);

1409:   KSPResetViewers(ksp);

1411:   ksp->setupstage = KSP_SETUP_NEW;
1412:   ksp->nmax = PETSC_DECIDE;
1413:   return(0);
1414: }

1416: /*@C
1417:    KSPDestroy - Destroys KSP context.

1419:    Collective on ksp

1421:    Input Parameter:
1422: .  ksp - iterative context obtained from KSPCreate()

1424:    Level: beginner

1426: .seealso: KSPCreate(), KSPSetUp(), KSPSolve(), KSP
1427: @*/
1428: PetscErrorCode  KSPDestroy(KSP *ksp)
1429: {
1431:   PC             pc;

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

1438:   PetscObjectSAWsViewOff((PetscObject)*ksp);

1440:   /*
1441:    Avoid a cascading call to PCReset(ksp->pc) from the following call:
1442:    PCReset() shouldn't be called from KSPDestroy() as it is unprotected by pc's
1443:    refcount (and may be shared, e.g., by other ksps).
1444:    */
1445:   pc         = (*ksp)->pc;
1446:   (*ksp)->pc = NULL;
1447:   KSPReset((*ksp));
1448:   (*ksp)->pc = pc;
1449:   if ((*ksp)->ops->destroy) {(*(*ksp)->ops->destroy)(*ksp);}

1451:   if ((*ksp)->transpose.use_explicittranspose) {
1452:     MatDestroy(&(*ksp)->transpose.AT);
1453:     MatDestroy(&(*ksp)->transpose.BT);
1454:     (*ksp)->transpose.reuse_transpose = PETSC_FALSE;
1455:   }

1457:   KSPGuessDestroy(&(*ksp)->guess);
1458:   DMDestroy(&(*ksp)->dm);
1459:   PCDestroy(&(*ksp)->pc);
1460:   PetscFree((*ksp)->res_hist_alloc);
1461:   PetscFree((*ksp)->err_hist_alloc);
1462:   if ((*ksp)->convergeddestroy) {
1463:     (*(*ksp)->convergeddestroy)((*ksp)->cnvP);
1464:   }
1465:   KSPMonitorCancel((*ksp));
1466:   KSPConvergedReasonViewCancel((*ksp));
1467:   PetscViewerDestroy(&(*ksp)->eigviewer);
1468:   PetscHeaderDestroy(ksp);
1469:   return(0);
1470: }

1472: /*@
1473:     KSPSetPCSide - Sets the preconditioning side.

1475:     Logically Collective on ksp

1477:     Input Parameter:
1478: .   ksp - iterative context obtained from KSPCreate()

1480:     Output Parameter:
1481: .   side - the preconditioning side, where side is one of
1482: .vb
1483:       PC_LEFT - left preconditioning (default)
1484:       PC_RIGHT - right preconditioning
1485:       PC_SYMMETRIC - symmetric preconditioning
1486: .ve

1488:     Options Database Keys:
1489: .   -ksp_pc_side <right,left,symmetric>

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

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

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

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

1502:     Level: intermediate

1504: .seealso: KSPGetPCSide(), KSPSetNormType(), KSPGetNormType(), KSP
1505: @*/
1506: PetscErrorCode  KSPSetPCSide(KSP ksp,PCSide side)
1507: {
1511:   ksp->pc_side = ksp->pc_side_set = side;
1512:   return(0);
1513: }

1515: /*@
1516:     KSPGetPCSide - Gets the preconditioning side.

1518:     Not Collective

1520:     Input Parameter:
1521: .   ksp - iterative context obtained from KSPCreate()

1523:     Output Parameter:
1524: .   side - the preconditioning side, where side is one of
1525: .vb
1526:       PC_LEFT - left preconditioning (default)
1527:       PC_RIGHT - right preconditioning
1528:       PC_SYMMETRIC - symmetric preconditioning
1529: .ve

1531:     Level: intermediate

1533: .seealso: KSPSetPCSide(), KSP
1534: @*/
1535: PetscErrorCode  KSPGetPCSide(KSP ksp,PCSide *side)
1536: {

1542:   KSPSetUpNorms_Private(ksp,PETSC_TRUE,&ksp->normtype,&ksp->pc_side);
1543:   *side = ksp->pc_side;
1544:   return(0);
1545: }

1547: /*@
1548:    KSPGetTolerances - Gets the relative, absolute, divergence, and maximum
1549:    iteration tolerances used by the default KSP convergence tests.

1551:    Not Collective

1553:    Input Parameter:
1554: .  ksp - the Krylov subspace context

1556:    Output Parameters:
1557: +  rtol - the relative convergence tolerance
1558: .  abstol - the absolute convergence tolerance
1559: .  dtol - the divergence tolerance
1560: -  maxits - maximum number of iterations

1562:    Notes:
1563:    The user can specify NULL for any parameter that is not needed.

1565:    Level: intermediate

1567:            maximum, iterations

1569: .seealso: KSPSetTolerances(), KSP
1570: @*/
1571: PetscErrorCode  KSPGetTolerances(KSP ksp,PetscReal *rtol,PetscReal *abstol,PetscReal *dtol,PetscInt *maxits)
1572: {
1575:   if (abstol) *abstol = ksp->abstol;
1576:   if (rtol) *rtol = ksp->rtol;
1577:   if (dtol) *dtol = ksp->divtol;
1578:   if (maxits) *maxits = ksp->max_it;
1579:   return(0);
1580: }

1582: /*@
1583:    KSPSetTolerances - Sets the relative, absolute, divergence, and maximum
1584:    iteration tolerances used by the default KSP convergence testers.

1586:    Logically Collective on ksp

1588:    Input Parameters:
1589: +  ksp - the Krylov subspace context
1590: .  rtol - the relative convergence tolerance, relative decrease in the (possibly preconditioned) residual norm
1591: .  abstol - the absolute convergence tolerance   absolute size of the (possibly preconditioned) residual norm
1592: .  dtol - the divergence tolerance,   amount (possibly preconditioned) residual norm can increase before KSPConvergedDefault() concludes that the method is diverging
1593: -  maxits - maximum number of iterations to use

1595:    Options Database Keys:
1596: +  -ksp_atol <abstol> - Sets abstol
1597: .  -ksp_rtol <rtol> - Sets rtol
1598: .  -ksp_divtol <dtol> - Sets dtol
1599: -  -ksp_max_it <maxits> - Sets maxits

1601:    Notes:
1602:    Use PETSC_DEFAULT to retain the default value of any of the tolerances.

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

1607:    Level: intermediate

1609:            convergence, maximum, iterations

1611: .seealso: KSPGetTolerances(), KSPConvergedDefault(), KSPSetConvergenceTest(), KSP
1612: @*/
1613: PetscErrorCode  KSPSetTolerances(KSP ksp,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt maxits)
1614: {

1622:   if (rtol != PETSC_DEFAULT) {
1623:     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);
1624:     ksp->rtol = rtol;
1625:   }
1626:   if (abstol != PETSC_DEFAULT) {
1627:     if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
1628:     ksp->abstol = abstol;
1629:   }
1630:   if (dtol != PETSC_DEFAULT) {
1631:     if (dtol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Divergence tolerance %g must be larger than 1.0",(double)dtol);
1632:     ksp->divtol = dtol;
1633:   }
1634:   if (maxits != PETSC_DEFAULT) {
1635:     if (maxits < 0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxits);
1636:     ksp->max_it = maxits;
1637:   }
1638:   return(0);
1639: }

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

1646:    Logically Collective on ksp

1648:    Input Parameters:
1649: +  ksp - iterative context obtained from KSPCreate()
1650: -  flg - PETSC_TRUE indicates the guess is non-zero, PETSC_FALSE indicates the guess is zero

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

1655:    Level: beginner

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

1660: .seealso: KSPGetInitialGuessNonzero(), KSPSetGuessType(), KSPGuessType, KSP
1661: @*/
1662: PetscErrorCode  KSPSetInitialGuessNonzero(KSP ksp,PetscBool flg)
1663: {
1667:   ksp->guess_zero = (PetscBool) !(int)flg;
1668:   return(0);
1669: }

1671: /*@
1672:    KSPGetInitialGuessNonzero - Determines whether the KSP solver is using
1673:    a zero initial guess.

1675:    Not Collective

1677:    Input Parameter:
1678: .  ksp - iterative context obtained from KSPCreate()

1680:    Output Parameter:
1681: .  flag - PETSC_TRUE if guess is nonzero, else PETSC_FALSE

1683:    Level: intermediate

1685: .seealso: KSPSetInitialGuessNonzero(), KSP
1686: @*/
1687: PetscErrorCode  KSPGetInitialGuessNonzero(KSP ksp,PetscBool  *flag)
1688: {
1692:   if (ksp->guess_zero) *flag = PETSC_FALSE;
1693:   else *flag = PETSC_TRUE;
1694:   return(0);
1695: }

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

1700:    Logically Collective on ksp

1702:    Input Parameters:
1703: +  ksp - iterative context obtained from KSPCreate()
1704: -  flg - PETSC_TRUE indicates you want the error generated

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

1709:    Level: intermediate

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


1716: .seealso: KSPGetErrorIfNotConverged(), KSP
1717: @*/
1718: PetscErrorCode  KSPSetErrorIfNotConverged(KSP ksp,PetscBool flg)
1719: {
1723:   ksp->errorifnotconverged = flg;
1724:   return(0);
1725: }

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

1730:    Not Collective

1732:    Input Parameter:
1733: .  ksp - iterative context obtained from KSPCreate()

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

1738:    Level: intermediate

1740: .seealso: KSPSetErrorIfNotConverged(), KSP
1741: @*/
1742: PetscErrorCode  KSPGetErrorIfNotConverged(KSP ksp,PetscBool  *flag)
1743: {
1747:   *flag = ksp->errorifnotconverged;
1748:   return(0);
1749: }

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

1754:    Logically Collective on ksp

1756:    Input Parameters:
1757: +  ksp - iterative context obtained from KSPCreate()
1758: -  flg - PETSC_TRUE or PETSC_FALSE

1760:    Level: advanced

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

1764: .seealso: KSPGetInitialGuessKnoll(), KSPSetInitialGuessNonzero(), KSPGetInitialGuessNonzero(), KSP
1765: @*/
1766: PetscErrorCode  KSPSetInitialGuessKnoll(KSP ksp,PetscBool flg)
1767: {
1771:   ksp->guess_knoll = flg;
1772:   return(0);
1773: }

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

1779:    Not Collective

1781:    Input Parameter:
1782: .  ksp - iterative context obtained from KSPCreate()

1784:    Output Parameter:
1785: .  flag - PETSC_TRUE if using Knoll trick, else PETSC_FALSE

1787:    Level: advanced

1789: .seealso: KSPSetInitialGuessKnoll(), KSPSetInitialGuessNonzero(), KSPGetInitialGuessNonzero(), KSP
1790: @*/
1791: PetscErrorCode  KSPGetInitialGuessKnoll(KSP ksp,PetscBool  *flag)
1792: {
1796:   *flag = ksp->guess_knoll;
1797:   return(0);
1798: }

1800: /*@
1801:    KSPGetComputeSingularValues - Gets the flag indicating whether the extreme singular
1802:    values will be calculated via a Lanczos or Arnoldi process as the linear
1803:    system is solved.

1805:    Not Collective

1807:    Input Parameter:
1808: .  ksp - iterative context obtained from KSPCreate()

1810:    Output Parameter:
1811: .  flg - PETSC_TRUE or PETSC_FALSE

1813:    Options Database Key:
1814: .  -ksp_monitor_singular_value - Activates KSPSetComputeSingularValues()

1816:    Notes:
1817:    Currently this option is not valid for all iterative methods.

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

1823:    Level: advanced

1825: .seealso: KSPComputeExtremeSingularValues(), KSPMonitorSingularValue(), KSP
1826: @*/
1827: PetscErrorCode  KSPGetComputeSingularValues(KSP ksp,PetscBool  *flg)
1828: {
1832:   *flg = ksp->calc_sings;
1833:   return(0);
1834: }

1836: /*@
1837:    KSPSetComputeSingularValues - Sets a flag so that the extreme singular
1838:    values will be calculated via a Lanczos or Arnoldi process as the linear
1839:    system is solved.

1841:    Logically Collective on ksp

1843:    Input Parameters:
1844: +  ksp - iterative context obtained from KSPCreate()
1845: -  flg - PETSC_TRUE or PETSC_FALSE

1847:    Options Database Key:
1848: .  -ksp_monitor_singular_value - Activates KSPSetComputeSingularValues()

1850:    Notes:
1851:    Currently this option is not valid for all iterative methods.

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

1857:    Level: advanced

1859: .seealso: KSPComputeExtremeSingularValues(), KSPMonitorSingularValue(), KSP
1860: @*/
1861: PetscErrorCode  KSPSetComputeSingularValues(KSP ksp,PetscBool flg)
1862: {
1866:   ksp->calc_sings = flg;
1867:   return(0);
1868: }

1870: /*@
1871:    KSPGetComputeEigenvalues - Gets the flag indicating that the extreme eigenvalues
1872:    values will be calculated via a Lanczos or Arnoldi process as the linear
1873:    system is solved.

1875:    Not Collective

1877:    Input Parameter:
1878: .  ksp - iterative context obtained from KSPCreate()

1880:    Output Parameter:
1881: .  flg - PETSC_TRUE or PETSC_FALSE

1883:    Notes:
1884:    Currently this option is not valid for all iterative methods.

1886:    Level: advanced

1888: .seealso: KSPComputeEigenvalues(), KSPComputeEigenvaluesExplicitly(), KSP
1889: @*/
1890: PetscErrorCode  KSPGetComputeEigenvalues(KSP ksp,PetscBool  *flg)
1891: {
1895:   *flg = ksp->calc_sings;
1896:   return(0);
1897: }

1899: /*@
1900:    KSPSetComputeEigenvalues - Sets a flag so that the extreme eigenvalues
1901:    values will be calculated via a Lanczos or Arnoldi process as the linear
1902:    system is solved.

1904:    Logically Collective on ksp

1906:    Input Parameters:
1907: +  ksp - iterative context obtained from KSPCreate()
1908: -  flg - PETSC_TRUE or PETSC_FALSE

1910:    Notes:
1911:    Currently this option is not valid for all iterative methods.

1913:    Level: advanced

1915: .seealso: KSPComputeEigenvalues(), KSPComputeEigenvaluesExplicitly(), KSP
1916: @*/
1917: PetscErrorCode  KSPSetComputeEigenvalues(KSP ksp,PetscBool flg)
1918: {
1922:   ksp->calc_sings = flg;
1923:   return(0);
1924: }

1926: /*@
1927:    KSPSetComputeRitz - Sets a flag so that the Ritz or harmonic Ritz pairs
1928:    will be calculated via a Lanczos or Arnoldi process as the linear
1929:    system is solved.

1931:    Logically Collective on ksp

1933:    Input Parameters:
1934: +  ksp - iterative context obtained from KSPCreate()
1935: -  flg - PETSC_TRUE or PETSC_FALSE

1937:    Notes:
1938:    Currently this option is only valid for the GMRES method.

1940:    Level: advanced

1942: .seealso: KSPComputeRitz(), KSP
1943: @*/
1944: PetscErrorCode  KSPSetComputeRitz(KSP ksp, PetscBool flg)
1945: {
1949:   ksp->calc_ritz = flg;
1950:   return(0);
1951: }

1953: /*@
1954:    KSPGetRhs - Gets the right-hand-side vector for the linear system to
1955:    be solved.

1957:    Not Collective

1959:    Input Parameter:
1960: .  ksp - iterative context obtained from KSPCreate()

1962:    Output Parameter:
1963: .  r - right-hand-side vector

1965:    Level: developer

1967: .seealso: KSPGetSolution(), KSPSolve(), KSP
1968: @*/
1969: PetscErrorCode  KSPGetRhs(KSP ksp,Vec *r)
1970: {
1974:   *r = ksp->vec_rhs;
1975:   return(0);
1976: }

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

1983:    Not Collective

1985:    Input Parameters:
1986: .  ksp - iterative context obtained from KSPCreate()

1988:    Output Parameters:
1989: .  v - solution vector

1991:    Level: developer

1993: .seealso: KSPGetRhs(),  KSPBuildSolution(), KSPSolve(), KSP
1994: @*/
1995: PetscErrorCode  KSPGetSolution(KSP ksp,Vec *v)
1996: {
2000:   *v = ksp->vec_sol;
2001:   return(0);
2002: }

2004: /*@
2005:    KSPSetPC - Sets the preconditioner to be used to calculate the
2006:    application of the preconditioner on a vector.

2008:    Collective on ksp

2010:    Input Parameters:
2011: +  ksp - iterative context obtained from KSPCreate()
2012: -  pc   - the preconditioner object (can be NULL)

2014:    Notes:
2015:    Use KSPGetPC() to retrieve the preconditioner context.

2017:    Level: developer

2019: .seealso: KSPGetPC(), KSP
2020: @*/
2021: PetscErrorCode  KSPSetPC(KSP ksp,PC pc)
2022: {

2027:   if (pc) {
2030:   }
2031:   PetscObjectReference((PetscObject)pc);
2032:   PCDestroy(&ksp->pc);
2033:   ksp->pc = pc;
2034:   PetscLogObjectParent((PetscObject)ksp,(PetscObject)ksp->pc);
2035:   return(0);
2036: }

2038: /*@
2039:    KSPGetPC - Returns a pointer to the preconditioner context
2040:    set with KSPSetPC().

2042:    Not Collective

2044:    Input Parameters:
2045: .  ksp - iterative context obtained from KSPCreate()

2047:    Output Parameter:
2048: .  pc - preconditioner context

2050:    Level: developer

2052: .seealso: KSPSetPC(), KSP
2053: @*/
2054: PetscErrorCode  KSPGetPC(KSP ksp,PC *pc)
2055: {

2061:   if (!ksp->pc) {
2062:     PCCreate(PetscObjectComm((PetscObject)ksp),&ksp->pc);
2063:     PetscObjectIncrementTabLevel((PetscObject)ksp->pc,(PetscObject)ksp,0);
2064:     PetscLogObjectParent((PetscObject)ksp,(PetscObject)ksp->pc);
2065:     PetscObjectSetOptions((PetscObject)ksp->pc,((PetscObject)ksp)->options);
2066:   }
2067:   *pc = ksp->pc;
2068:   return(0);
2069: }

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

2074:    Collective on ksp

2076:    Input Parameters:
2077: +  ksp - iterative context obtained from KSPCreate()
2078: .  it - iteration number
2079: -  rnorm - relative norm of the residual

2081:    Notes:
2082:    This routine is called by the KSP implementations.
2083:    It does not typically need to be called by the user.

2085:    Level: developer

2087: .seealso: KSPMonitorSet()
2088: @*/
2089: PetscErrorCode KSPMonitor(KSP ksp,PetscInt it,PetscReal rnorm)
2090: {
2091:   PetscInt       i, n = ksp->numbermonitors;

2095:   for (i=0; i<n; i++) {
2096:     (*ksp->monitor[i])(ksp,it,rnorm,ksp->monitorcontext[i]);
2097:   }
2098:   return(0);
2099: }

2101: /*@C
2102:    KSPMonitorSet - Sets an ADDITIONAL function to be called at every iteration to monitor
2103:    the residual/error etc.

2105:    Logically Collective on ksp

2107:    Input Parameters:
2108: +  ksp - iterative context obtained from KSPCreate()
2109: .  monitor - pointer to function (if this is NULL, it turns off monitoring
2110: .  mctx    - [optional] context for private data for the
2111:              monitor routine (use NULL if no context is desired)
2112: -  monitordestroy - [optional] routine that frees monitor context
2113:           (may be NULL)

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

2118: +  ksp - iterative context obtained from KSPCreate()
2119: .  it - iteration number
2120: .  rnorm - (estimated) 2-norm of (preconditioned) residual
2121: -  mctx  - optional monitoring context, as set by KSPMonitorSet()

2123:    Options Database Keys:
2124: +    -ksp_monitor               - sets KSPMonitorResidual()
2125: .    -ksp_monitor draw          - sets KSPMonitorResidualDraw() and plots residual
2126: .    -ksp_monitor draw::draw_lg - sets KSPMonitorResidualDrawLG() and plots residual
2127: .    -ksp_monitor_pause_final   - Pauses any graphics when the solve finishes (only works for internal monitors)
2128: .    -ksp_monitor_true_residual - sets KSPMonitorTrueResidual()
2129: .    -ksp_monitor_true_residual draw::draw_lg - sets KSPMonitorTrueResidualDrawLG() and plots residual
2130: .    -ksp_monitor_max           - sets KSPMonitorTrueResidualMax()
2131: .    -ksp_monitor_singular_value - sets KSPMonitorSingularValue()
2132: -    -ksp_monitor_cancel - cancels all monitors that have
2133:                           been hardwired into a code by
2134:                           calls to KSPMonitorSet(), but
2135:                           does not cancel those set via
2136:                           the options database.

2138:    Notes:
2139:    The default is to do nothing.  To print the residual, or preconditioned
2140:    residual if KSPSetNormType(ksp,KSP_NORM_PRECONDITIONED) was called, use
2141:    KSPMonitorResidual() as the monitoring routine, with a ASCII viewer as the
2142:    context.

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

2148:    Fortran Notes:
2149:     Only a single monitor function can be set for each KSP object

2151:    Level: beginner

2153: .seealso: KSPMonitorResidual(), KSPMonitorCancel(), KSP
2154: @*/
2155: PetscErrorCode  KSPMonitorSet(KSP ksp,PetscErrorCode (*monitor)(KSP,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
2156: {
2157:   PetscInt       i;
2159:   PetscBool      identical;

2163:   for (i=0; i<ksp->numbermonitors;i++) {
2164:     PetscMonitorCompare((PetscErrorCode (*)(void))monitor,mctx,monitordestroy,(PetscErrorCode (*)(void))ksp->monitor[i],ksp->monitorcontext[i],ksp->monitordestroy[i],&identical);
2165:     if (identical) return(0);
2166:   }
2167:   if (ksp->numbermonitors >= MAXKSPMONITORS) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Too many KSP monitors set");
2168:   ksp->monitor[ksp->numbermonitors]          = monitor;
2169:   ksp->monitordestroy[ksp->numbermonitors]   = monitordestroy;
2170:   ksp->monitorcontext[ksp->numbermonitors++] = (void*)mctx;
2171:   return(0);
2172: }

2174: /*@
2175:    KSPMonitorCancel - Clears all monitors for a KSP object.

2177:    Logically Collective on ksp

2179:    Input Parameters:
2180: .  ksp - iterative context obtained from KSPCreate()

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

2187:    Level: intermediate

2189: .seealso: KSPMonitorResidual(), KSPMonitorSet(), KSP
2190: @*/
2191: PetscErrorCode  KSPMonitorCancel(KSP ksp)
2192: {
2194:   PetscInt       i;

2198:   for (i=0; i<ksp->numbermonitors; i++) {
2199:     if (ksp->monitordestroy[i]) {
2200:       (*ksp->monitordestroy[i])(&ksp->monitorcontext[i]);
2201:     }
2202:   }
2203:   ksp->numbermonitors = 0;
2204:   return(0);
2205: }

2207: /*@C
2208:    KSPGetMonitorContext - Gets the monitoring context, as set by
2209:    KSPMonitorSet() for the FIRST monitor only.

2211:    Not Collective

2213:    Input Parameter:
2214: .  ksp - iterative context obtained from KSPCreate()

2216:    Output Parameter:
2217: .  ctx - monitoring context

2219:    Level: intermediate

2221: .seealso: KSPMonitorResidual(), KSP
2222: @*/
2223: PetscErrorCode  KSPGetMonitorContext(KSP ksp,void **ctx)
2224: {
2227:   *ctx =      (ksp->monitorcontext[0]);
2228:   return(0);
2229: }

2231: /*@
2232:    KSPSetResidualHistory - Sets the array used to hold the residual history.
2233:    If set, this array will contain the residual norms computed at each
2234:    iteration of the solver.

2236:    Not Collective

2238:    Input Parameters:
2239: +  ksp - iterative context obtained from KSPCreate()
2240: .  a   - array to hold history
2241: .  na  - size of a
2242: -  reset - PETSC_TRUE indicates the history counter is reset to zero
2243:            for each new linear solve

2245:    Level: advanced

2247:    Notes:
2248:    If provided, he array is NOT freed by PETSc so the user needs to keep track of it and destroy once the KSP object is destroyed.
2249:    If 'a' is NULL then space is allocated for the history. If 'na' PETSC_DECIDE or PETSC_DEFAULT then a
2250:    default array of length 10000 is allocated.

2252: .seealso: KSPGetResidualHistory(), KSP

2254: @*/
2255: PetscErrorCode KSPSetResidualHistory(KSP ksp,PetscReal a[],PetscInt na,PetscBool reset)
2256: {


2262:   PetscFree(ksp->res_hist_alloc);
2263:   if (na != PETSC_DECIDE && na != PETSC_DEFAULT && a) {
2264:     ksp->res_hist     = a;
2265:     ksp->res_hist_max = na;
2266:   } else {
2267:     if (na != PETSC_DECIDE && na != PETSC_DEFAULT) ksp->res_hist_max = na;
2268:     else                                           ksp->res_hist_max = 10000; /* like default ksp->max_it */
2269:     PetscCalloc1(ksp->res_hist_max,&ksp->res_hist_alloc);

2271:     ksp->res_hist = ksp->res_hist_alloc;
2272:   }
2273:   ksp->res_hist_len   = 0;
2274:   ksp->res_hist_reset = reset;
2275:   return(0);
2276: }

2278: /*@C
2279:    KSPGetResidualHistory - Gets the array used to hold the residual history
2280:    and the number of residuals it contains.

2282:    Not Collective

2284:    Input Parameter:
2285: .  ksp - iterative context obtained from KSPCreate()

2287:    Output Parameters:
2288: +  a   - pointer to array to hold history (or NULL)
2289: -  na  - number of used entries in a (or NULL)

2291:    Level: advanced

2293:    Notes:
2294:      This array is borrowed and should not be freed by the caller.
2295:      Can only be called after a KSPSetResidualHistory() otherwise a and na are set to zero

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

2303: .seealso: KSPSetResidualHistory(), KSP

2305: @*/
2306: PetscErrorCode KSPGetResidualHistory(KSP ksp, const PetscReal *a[],PetscInt *na)
2307: {
2310:   if (a) *a = ksp->res_hist;
2311:   if (na) *na = ksp->res_hist_len;
2312:   return(0);
2313: }

2315: /*@
2316:   KSPSetErrorHistory - Sets the array used to hold the error history. If set, this array will contain the error norms computed at each iteration of the solver.

2318:   Not Collective

2320:   Input Parameters:
2321: + ksp   - iterative context obtained from KSPCreate()
2322: . a     - array to hold history
2323: . na    - size of a
2324: - reset - PETSC_TRUE indicates the history counter is reset to zero for each new linear solve

2326:   Level: advanced

2328:   Notes:
2329:   If provided, the array is NOT freed by PETSc so the user needs to keep track of it and destroy once the KSP object is destroyed.
2330:   If 'a' is NULL then space is allocated for the history. If 'na' PETSC_DECIDE or PETSC_DEFAULT then a default array of length 10000 is allocated.

2332: .seealso: KSPGetErrorHistory(), KSPSetResidualHistory(), KSP
2333: @*/
2334: PetscErrorCode KSPSetErrorHistory(KSP ksp, PetscReal a[], PetscInt na, PetscBool reset)
2335: {


2341:   PetscFree(ksp->err_hist_alloc);
2342:   if (na != PETSC_DECIDE && na != PETSC_DEFAULT && a) {
2343:     ksp->err_hist     = a;
2344:     ksp->err_hist_max = na;
2345:   } else {
2346:     if (na != PETSC_DECIDE && na != PETSC_DEFAULT) ksp->err_hist_max = na;
2347:     else                                           ksp->err_hist_max = 10000; /* like default ksp->max_it */
2348:     PetscCalloc1(ksp->err_hist_max, &ksp->err_hist_alloc);

2350:     ksp->err_hist = ksp->err_hist_alloc;
2351:   }
2352:   ksp->err_hist_len   = 0;
2353:   ksp->err_hist_reset = reset;
2354:   return(0);
2355: }

2357: /*@C
2358:   KSPGetErrorHistory - Gets the array used to hold the error history and the number of residuals it contains.

2360:   Not Collective

2362:   Input Parameter:
2363: . ksp - iterative context obtained from KSPCreate()

2365:   Output Parameters:
2366: + a  - pointer to array to hold history (or NULL)
2367: - na - number of used entries in a (or NULL)

2369:   Level: advanced

2371:   Notes:
2372:   This array is borrowed and should not be freed by the caller.
2373:   Can only be called after a KSPSetErrorHistory() otherwise a and na are set to zero
2374:   The Fortran version of this routine has a calling sequence
2375: $   call KSPGetErrorHistory(KSP ksp, integer na, integer ierr)
2376:   note that you have passed a Fortran array into KSPSetErrorHistory() and you need
2377:   to access the residual values from this Fortran array you provided. Only the na (number of
2378:   residual norms currently held) is set.

2380: .seealso: KSPSetErrorHistory(), KSPGetResidualHistory(), KSP
2381: @*/
2382: PetscErrorCode KSPGetErrorHistory(KSP ksp, const PetscReal *a[], PetscInt *na)
2383: {
2386:   if (a)  *a  = ksp->err_hist;
2387:   if (na) *na = ksp->err_hist_len;
2388:   return(0);
2389: }

2391: /*
2392:   KSPComputeConvergenceRate - Compute the convergence rate for the iteration

2394:   Not collective

2396:   Input Parameter:
2397: . ksp - The KSP

2399:   Output Parameters:
2400: + cr   - The residual contraction rate
2401: . rRsq - The coefficient of determination, R^2, indicating the linearity of the data
2402: . ce   - The error contraction rate
2403: - eRsq - The coefficient of determination, R^2, indicating the linearity of the data

2405:   Note:
2406:   Suppose that the residual is reduced linearly, $r_k = c^k r_0$, which means $log r_k = log r_0 + k log c$. After linear regression,
2407:   the slope is $\log c$. The coefficient of determination is given by $1 - \frac{\sum_i (y_i - f(x_i))^2}{\sum_i (y_i - \bar y)}$,
2408:   see also https://en.wikipedia.org/wiki/Coefficient_of_determination

2410:   Level: advanced

2412: .seealso: KSPConvergedRateView()
2413: */
2414: PetscErrorCode KSPComputeConvergenceRate(KSP ksp, PetscReal *cr, PetscReal *rRsq, PetscReal *ce, PetscReal *eRsq)
2415: {
2416:   PetscReal      const *hist;
2417:   PetscReal      *x, *y, slope, intercept, mean = 0.0, var = 0.0, res = 0.0;
2418:   PetscInt       n, k;

2422:   if (cr || rRsq) {
2423:     KSPGetResidualHistory(ksp, &hist, &n);
2424:     if (!n) {
2425:       if (cr)   *cr   =  0.0;
2426:       if (rRsq) *rRsq = -1.0;
2427:     } else {
2428:       PetscMalloc2(n, &x, n, &y);
2429:       for (k = 0; k < n; ++k) {
2430:         x[k] = k;
2431:         y[k] = PetscLogReal(hist[k]);
2432:         mean += y[k];
2433:       }
2434:       mean /= n;
2435:       PetscLinearRegression(n, x, y, &slope, &intercept);
2436:       for (k = 0; k < n; ++k) {
2437:         res += PetscSqr(y[k] - (slope*x[k] + intercept));
2438:         var += PetscSqr(y[k] - mean);
2439:       }
2440:       PetscFree2(x, y);
2441:       if (cr)   *cr   = PetscExpReal(slope);
2442:       if (rRsq) *rRsq = var < PETSC_MACHINE_EPSILON ? 0.0 : 1.0 - (res / var);
2443:     }
2444:   }
2445:   if (ce || eRsq) {
2446:     KSPGetErrorHistory(ksp, &hist, &n);
2447:     if (!n) {
2448:       if (ce)   *ce   =  0.0;
2449:       if (eRsq) *eRsq = -1.0;
2450:     } else {
2451:       PetscMalloc2(n, &x, n, &y);
2452:       for (k = 0; k < n; ++k) {
2453:         x[k] = k;
2454:         y[k] = PetscLogReal(hist[k]);
2455:         mean += y[k];
2456:       }
2457:       mean /= n;
2458:       PetscLinearRegression(n, x, y, &slope, &intercept);
2459:       for (k = 0; k < n; ++k) {
2460:         res += PetscSqr(y[k] - (slope*x[k] + intercept));
2461:         var += PetscSqr(y[k] - mean);
2462:       }
2463:       PetscFree2(x, y);
2464:       if (ce)   *ce   = PetscExpReal(slope);
2465:       if (eRsq) *eRsq = var < PETSC_MACHINE_EPSILON ? 0.0 : 1.0 - (res / var);
2466:     }
2467:   }
2468:   return(0);
2469: }

2471: /*@C
2472:    KSPSetConvergenceTest - Sets the function to be used to determine
2473:    convergence.

2475:    Logically Collective on ksp

2477:    Input Parameters:
2478: +  ksp - iterative context obtained from KSPCreate()
2479: .  converge - pointer to the function
2480: .  cctx    - context for private data for the convergence routine (may be null)
2481: -  destroy - a routine for destroying the context (may be null)

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

2486: +  ksp - iterative context obtained from KSPCreate()
2487: .  it - iteration number
2488: .  rnorm - (estimated) 2-norm of (preconditioned) residual
2489: .  reason - the reason why it has converged or diverged
2490: -  cctx  - optional convergence context, as set by KSPSetConvergenceTest()


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

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

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

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

2507:    Level: advanced

2509: .seealso: KSPConvergedDefault(), KSPGetConvergenceContext(), KSPSetTolerances(), KSP, KSPGetConvergenceTest(), KSPGetAndClearConvergenceTest()
2510: @*/
2511: PetscErrorCode  KSPSetConvergenceTest(KSP ksp,PetscErrorCode (*converge)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
2512: {

2517:   if (ksp->convergeddestroy) {
2518:     (*ksp->convergeddestroy)(ksp->cnvP);
2519:   }
2520:   ksp->converged        = converge;
2521:   ksp->convergeddestroy = destroy;
2522:   ksp->cnvP             = (void*)cctx;
2523:   return(0);
2524: }

2526: /*@C
2527:    KSPGetConvergenceTest - Gets the function to be used to determine
2528:    convergence.

2530:    Logically Collective on ksp

2532:    Input Parameter:
2533: .   ksp - iterative context obtained from KSPCreate()

2535:    Output Parameter:
2536: +  converge - pointer to convergence test function
2537: .  cctx    - context for private data for the convergence routine (may be null)
2538: -  destroy - a routine for destroying the context (may be null)

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

2543: +  ksp - iterative context obtained from KSPCreate()
2544: .  it - iteration number
2545: .  rnorm - (estimated) 2-norm of (preconditioned) residual
2546: .  reason - the reason why it has converged or diverged
2547: -  cctx  - optional convergence context, as set by KSPSetConvergenceTest()

2549:    Level: advanced

2551: .seealso: KSPConvergedDefault(), KSPGetConvergenceContext(), KSPSetTolerances(), KSP, KSPSetConvergenceTest(), KSPGetAndClearConvergenceTest()
2552: @*/
2553: PetscErrorCode  KSPGetConvergenceTest(KSP ksp,PetscErrorCode (**converge)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void **cctx,PetscErrorCode (**destroy)(void*))
2554: {
2557:   if (converge) *converge = ksp->converged;
2558:   if (destroy)  *destroy  = ksp->convergeddestroy;
2559:   if (cctx)     *cctx     = ksp->cnvP;
2560:   return(0);
2561: }

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

2566:    Logically Collective on ksp

2568:    Input Parameter:
2569: .   ksp - iterative context obtained from KSPCreate()

2571:    Output Parameter:
2572: +  converge - pointer to convergence test function
2573: .  cctx    - context for private data for the convergence routine
2574: -  destroy - a routine for destroying the context

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

2579: +  ksp - iterative context obtained from KSPCreate()
2580: .  it - iteration number
2581: .  rnorm - (estimated) 2-norm of (preconditioned) residual
2582: .  reason - the reason why it has converged or diverged
2583: -  cctx  - optional convergence context, as set by KSPSetConvergenceTest()

2585:    Level: advanced

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

2591: .seealso: KSPConvergedDefault(), KSPGetConvergenceContext(), KSPSetTolerances(), KSP, KSPSetConvergenceTest(), KSPGetConvergenceTest()
2592: @*/
2593: PetscErrorCode  KSPGetAndClearConvergenceTest(KSP ksp,PetscErrorCode (**converge)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void **cctx,PetscErrorCode (**destroy)(void*))
2594: {
2597:   *converge             = ksp->converged;
2598:   *destroy              = ksp->convergeddestroy;
2599:   *cctx                 = ksp->cnvP;
2600:   ksp->converged        = NULL;
2601:   ksp->cnvP             = NULL;
2602:   ksp->convergeddestroy = NULL;
2603:   return(0);
2604: }

2606: /*@C
2607:    KSPGetConvergenceContext - Gets the convergence context set with
2608:    KSPSetConvergenceTest().

2610:    Not Collective

2612:    Input Parameter:
2613: .  ksp - iterative context obtained from KSPCreate()

2615:    Output Parameter:
2616: .  ctx - monitoring context

2618:    Level: advanced

2620: .seealso: KSPConvergedDefault(), KSPSetConvergenceTest(), KSP
2621: @*/
2622: PetscErrorCode  KSPGetConvergenceContext(KSP ksp,void **ctx)
2623: {
2626:   *ctx = ksp->cnvP;
2627:   return(0);
2628: }

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

2634:    Collective on ksp

2636:    Input Parameter:
2637: .  ctx - iterative context obtained from KSPCreate()

2639:    Output Parameter:
2640:    Provide exactly one of
2641: +  v - location to stash solution.
2642: -  V - the solution is returned in this location. This vector is created
2643:        internally. This vector should NOT be destroyed by the user with
2644:        VecDestroy().

2646:    Notes:
2647:    This routine can be used in one of two ways
2648: .vb
2649:       KSPBuildSolution(ksp,NULL,&V);
2650:    or
2651:       KSPBuildSolution(ksp,v,NULL); or KSPBuildSolution(ksp,v,&v);
2652: .ve
2653:    In the first case an internal vector is allocated to store the solution
2654:    (the user cannot destroy this vector). In the second case the solution
2655:    is generated in the vector that the user provides. Note that for certain
2656:    methods, such as KSPCG, the second case requires a copy of the solution,
2657:    while in the first case the call is essentially free since it simply
2658:    returns the vector where the solution already is stored. For some methods
2659:    like GMRES this is a reasonably expensive operation and should only be
2660:    used in truly needed.

2662:    Level: advanced

2664: .seealso: KSPGetSolution(), KSPBuildResidual(), KSP
2665: @*/
2666: PetscErrorCode  KSPBuildSolution(KSP ksp,Vec v,Vec *V)
2667: {

2672:   if (!V && !v) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONG,"Must provide either v or V");
2673:   if (!V) V = &v;
2674:   (*ksp->ops->buildsolution)(ksp,v,V);
2675:   return(0);
2676: }

2678: /*@C
2679:    KSPBuildResidual - Builds the residual in a vector provided.

2681:    Collective on ksp

2683:    Input Parameter:
2684: .  ksp - iterative context obtained from KSPCreate()

2686:    Output Parameters:
2687: +  v - optional location to stash residual.  If v is not provided,
2688:        then a location is generated.
2689: .  t - work vector.  If not provided then one is generated.
2690: -  V - the residual

2692:    Notes:
2693:    Regardless of whether or not v is provided, the residual is
2694:    returned in V.

2696:    Level: advanced

2698: .seealso: KSPBuildSolution()
2699: @*/
2700: PetscErrorCode  KSPBuildResidual(KSP ksp,Vec t,Vec v,Vec *V)
2701: {
2703:   PetscBool      flag = PETSC_FALSE;
2704:   Vec            w    = v,tt = t;

2708:   if (!w) {
2709:     VecDuplicate(ksp->vec_rhs,&w);
2710:     PetscLogObjectParent((PetscObject)ksp,(PetscObject)w);
2711:   }
2712:   if (!tt) {
2713:     VecDuplicate(ksp->vec_sol,&tt); flag = PETSC_TRUE;
2714:     PetscLogObjectParent((PetscObject)ksp,(PetscObject)tt);
2715:   }
2716:   (*ksp->ops->buildresidual)(ksp,tt,w,V);
2717:   if (flag) {VecDestroy(&tt);}
2718:   return(0);
2719: }

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

2725:    Logically Collective on ksp

2727:    Input Parameter:
2728: +  ksp - the KSP context
2729: -  scale - PETSC_TRUE or PETSC_FALSE

2731:    Options Database Key:
2732: +   -ksp_diagonal_scale -
2733: -   -ksp_diagonal_scale_fix - scale the matrix back AFTER the solve


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

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

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

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

2751:    Level: intermediate

2753: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScaleFix(), KSP
2754: @*/
2755: PetscErrorCode  KSPSetDiagonalScale(KSP ksp,PetscBool scale)
2756: {
2760:   ksp->dscale = scale;
2761:   return(0);
2762: }

2764: /*@
2765:    KSPGetDiagonalScale - Checks if KSP solver scales the matrix and
2766:                           right hand side

2768:    Not Collective

2770:    Input Parameter:
2771: .  ksp - the KSP context

2773:    Output Parameter:
2774: .  scale - PETSC_TRUE or PETSC_FALSE

2776:    Notes:
2777:     BE CAREFUL with this routine: it actually scales the matrix and right
2778:     hand side that define the system. After the system is solved the matrix
2779:     and right hand side remain scaled  unless you use KSPSetDiagonalScaleFix()

2781:    Level: intermediate

2783: .seealso: KSPSetDiagonalScale(), KSPSetDiagonalScaleFix(), KSP
2784: @*/
2785: PetscErrorCode  KSPGetDiagonalScale(KSP ksp,PetscBool  *scale)
2786: {
2790:   *scale = ksp->dscale;
2791:   return(0);
2792: }

2794: /*@
2795:    KSPSetDiagonalScaleFix - Tells KSP to diagonally scale the system
2796:      back after solving.

2798:    Logically Collective on ksp

2800:    Input Parameter:
2801: +  ksp - the KSP context
2802: -  fix - PETSC_TRUE to scale back after the system solve, PETSC_FALSE to not
2803:          rescale (default)

2805:    Notes:
2806:      Must be called after KSPSetDiagonalScale()

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

2813:    Level: intermediate

2815: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScale(), KSPGetDiagonalScaleFix(), KSP
2816: @*/
2817: PetscErrorCode  KSPSetDiagonalScaleFix(KSP ksp,PetscBool fix)
2818: {
2822:   ksp->dscalefix = fix;
2823:   return(0);
2824: }

2826: /*@
2827:    KSPGetDiagonalScaleFix - Determines if KSP diagonally scales the system
2828:      back after solving.

2830:    Not Collective

2832:    Input Parameter:
2833: .  ksp - the KSP context

2835:    Output Parameter:
2836: .  fix - PETSC_TRUE to scale back after the system solve, PETSC_FALSE to not
2837:          rescale (default)

2839:    Notes:
2840:      Must be called after KSPSetDiagonalScale()

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

2847:    Level: intermediate

2849: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScale(), KSPSetDiagonalScaleFix(), KSP
2850: @*/
2851: PetscErrorCode  KSPGetDiagonalScaleFix(KSP ksp,PetscBool  *fix)
2852: {
2856:   *fix = ksp->dscalefix;
2857:   return(0);
2858: }

2860: /*@C
2861:    KSPSetComputeOperators - set routine to compute the linear operators

2863:    Logically Collective

2865:    Input Arguments:
2866: +  ksp - the KSP context
2867: .  func - function to compute the operators
2868: -  ctx - optional context

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

2873: +  ksp - the KSP context
2874: .  A - the linear operator
2875: .  B - preconditioning matrix
2876: -  ctx - optional user-provided context

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

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

2884:    Level: beginner

2886: .seealso: KSPSetOperators(), KSPSetComputeRHS(), DMKSPSetComputeOperators(), KSPSetComputeInitialGuess()
2887: @*/
2888: PetscErrorCode KSPSetComputeOperators(KSP ksp,PetscErrorCode (*func)(KSP,Mat,Mat,void*),void *ctx)
2889: {
2891:   DM             dm;

2895:   KSPGetDM(ksp,&dm);
2896:   DMKSPSetComputeOperators(dm,func,ctx);
2897:   if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX;
2898:   return(0);
2899: }

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

2904:    Logically Collective

2906:    Input Arguments:
2907: +  ksp - the KSP context
2908: .  func - function to compute the right hand side
2909: -  ctx - optional context

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

2914: +  ksp - the KSP context
2915: .  b - right hand side of linear system
2916: -  ctx - optional user-provided context

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

2921:    Level: beginner

2923: .seealso: KSPSolve(), DMKSPSetComputeRHS(), KSPSetComputeOperators()
2924: @*/
2925: PetscErrorCode KSPSetComputeRHS(KSP ksp,PetscErrorCode (*func)(KSP,Vec,void*),void *ctx)
2926: {
2928:   DM             dm;

2932:   KSPGetDM(ksp,&dm);
2933:   DMKSPSetComputeRHS(dm,func,ctx);
2934:   return(0);
2935: }

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

2940:    Logically Collective

2942:    Input Arguments:
2943: +  ksp - the KSP context
2944: .  func - function to compute the initial guess
2945: -  ctx - optional context

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

2950: +  ksp - the KSP context
2951: .  x - solution vector
2952: -  ctx - optional user-provided context

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

2957:    Level: beginner

2959: .seealso: KSPSolve(), KSPSetComputeRHS(), KSPSetComputeOperators(), DMKSPSetComputeInitialGuess()
2960: @*/
2961: PetscErrorCode KSPSetComputeInitialGuess(KSP ksp,PetscErrorCode (*func)(KSP,Vec,void*),void *ctx)
2962: {
2964:   DM             dm;

2968:   KSPGetDM(ksp,&dm);
2969:   DMKSPSetComputeInitialGuess(dm,func,ctx);
2970:   return(0);
2971: }

2973: /*@
2974:    KSPSetUseExplicitTranspose - Determines if transpose the system explicitly
2975:    in KSPSolveTranspose.

2977:    Logically Collective on ksp

2979:    Input Parameter:
2980: .  ksp - the KSP context

2982:    Output Parameter:
2983: .  flg - PETSC_TRUE to transpose the system in KSPSolveTranspose, PETSC_FALSE to not
2984:          transpose (default)

2986:    Level: advanced

2988: .seealso: KSPSolveTranspose(), KSP
2989: @*/
2990: PetscErrorCode KSPSetUseExplicitTranspose(KSP ksp,PetscBool flg)
2991: {
2994:   ksp->transpose.use_explicittranspose = flg;
2995:   return(0);
2996: }