Actual source code: itfunc.c

petsc-3.5.2 2014-09-08
Report Typos and Errors
  2: /*
  3:       Interface KSP routines that the user calls.
  4: */

  6: #include <petsc-private/kspimpl.h>   /*I "petscksp.h" I*/
  7: #include <petscdm.h>

 11: /*@
 12:    KSPComputeExtremeSingularValues - Computes the extreme singular values
 13:    for the preconditioned operator. Called after or during KSPSolve().

 15:    Not Collective

 17:    Input Parameter:
 18: .  ksp - iterative context obtained from KSPCreate()

 20:    Output Parameters:
 21: .  emin, emax - extreme singular values

 23:    Options Database Keys:
 24: .  -ksp_compute_singularvalues - compute extreme singular values and print when KSPSolve completes.

 26:    Notes:
 27:    One must call KSPSetComputeSingularValues() before calling KSPSetUp()
 28:    (or use the option -ksp_compute_eigenvalues) in order for this routine to work correctly.

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

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

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

 41:    Level: advanced

 43: .keywords: KSP, compute, extreme, singular, values

 45: .seealso: KSPSetComputeSingularValues(), KSPMonitorSingularValue(), KSPComputeEigenvalues()
 46: @*/
 47: PetscErrorCode  KSPComputeExtremeSingularValues(KSP ksp,PetscReal *emax,PetscReal *emin)
 48: {

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

 57:   if (ksp->ops->computeextremesingularvalues) {
 58:     (*ksp->ops->computeextremesingularvalues)(ksp,emax,emin);
 59:   } else {
 60:     *emin = -1.0;
 61:     *emax = -1.0;
 62:   }
 63:   return(0);
 64: }

 68: /*@
 69:    KSPComputeEigenvalues - Computes the extreme eigenvalues for the
 70:    preconditioned operator. Called after or during KSPSolve().

 72:    Not Collective

 74:    Input Parameter:
 75: +  ksp - iterative context obtained from KSPCreate()
 76: -  n - size of arrays r and c. The number of eigenvalues computed (neig) will, in
 77:        general, be less than this.

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

 84:    Options Database Keys:
 85: +  -ksp_compute_eigenvalues - Prints eigenvalues to stdout
 86: -  -ksp_plot_eigenvalues - Plots eigenvalues in an x-window display

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

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

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

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

107:    Level: advanced

109: .keywords: KSP, compute, extreme, singular, values

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

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

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

134: /*@
135:    KSPSetUpOnBlocks - Sets up the preconditioner for each block in
136:    the block Jacobi, block Gauss-Seidel, and overlapping Schwarz
137:    methods.

139:    Collective on KSP

141:    Input Parameter:
142: .  ksp - the KSP context

144:    Notes:
145:    KSPSetUpOnBlocks() is a routine that the user can optinally call for
146:    more precise profiling (via -log_summary) of the setup phase for these
147:    block preconditioners.  If the user does not call KSPSetUpOnBlocks(),
148:    it will automatically be called from within KSPSolve().

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

153:    Level: advanced

155: .keywords: KSP, setup, blocks

157: .seealso: PCSetUpOnBlocks(), KSPSetUp(), PCSetUp()
158: @*/
159: PetscErrorCode  KSPSetUpOnBlocks(KSP ksp)
160: {

165:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
166:   PCSetUpOnBlocks(ksp->pc);
167:   return(0);
168: }

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

175:    Collective on KSP

177:    Input Parameters:
178: +  ksp   - iterative context obtained from KSPCreate()
179: -  flag - PETSC_TRUE to reuse the current preconditioner

181:    Level: intermediate

183: .keywords: KSP, setup

185: .seealso: KSPCreate(), KSPSolve(), KSPDestroy(), PCSetReusePreconditioner()
186: @*/
187: PetscErrorCode  KSPSetReusePreconditioner(KSP ksp,PetscBool flag)
188: {

193:   PCSetReusePreconditioner(ksp->pc,flag);
194:   return(0);
195: }

199: /*@
200:    KSPSetUp - Sets up the internal data structures for the
201:    later use of an iterative solver.

203:    Collective on KSP

205:    Input Parameter:
206: .  ksp   - iterative context obtained from KSPCreate()

208:    Level: developer

210: .keywords: KSP, setup

212: .seealso: KSPCreate(), KSPSolve(), KSPDestroy()
213: @*/
214: PetscErrorCode  KSPSetUp(KSP ksp)
215: {
217:   Mat            A,B;


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

225:   if (!((PetscObject)ksp)->type_name) {
226:     KSPSetType(ksp,KSPGMRES);
227:   }
228:   KSPSetUpNorms_Private(ksp,&ksp->normtype,&ksp->pc_side);

230:   if (ksp->dmActive && !ksp->setupstage) {
231:     /* first time in so build matrix and vector data structures using DM */
232:     if (!ksp->vec_rhs) {DMCreateGlobalVector(ksp->dm,&ksp->vec_rhs);}
233:     if (!ksp->vec_sol) {DMCreateGlobalVector(ksp->dm,&ksp->vec_sol);}
234:     DMCreateMatrix(ksp->dm,&A);
235:     KSPSetOperators(ksp,A,A);
236:     PetscObjectDereference((PetscObject)A);
237:   }

239:   if (ksp->dmActive) {
240:     DMKSP kdm;
241:     DMGetDMKSP(ksp->dm,&kdm);

243:     if (kdm->ops->computeinitialguess && ksp->setupstage != KSP_SETUP_NEWRHS) {
244:       /* only computes initial guess the first time through */
245:       (*kdm->ops->computeinitialguess)(ksp,ksp->vec_sol,kdm->initialguessctx);
246:       KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
247:     }
248:     if (kdm->ops->computerhs) {
249:       (*kdm->ops->computerhs)(ksp,ksp->vec_rhs,kdm->rhsctx);
250:     }

252:     if (ksp->setupstage != KSP_SETUP_NEWRHS) {
253:       if (kdm->ops->computeoperators) {
254:         KSPGetOperators(ksp,&A,&B);
255:         (*kdm->ops->computeoperators)(ksp,A,B,kdm->operatorsctx);
256:         KSPSetOperators(ksp,A,B);
257:       } else SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONGSTATE,"You called KSPSetDM() but did not use DMKSPSetComputeOperators() or KSPSetDMActive(dm,PETSC_FALSE);");
258:     }
259:   }

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

264:   switch (ksp->setupstage) {
265:   case KSP_SETUP_NEW:
266:     (*ksp->ops->setup)(ksp);
267:     break;
268:   case KSP_SETUP_NEWMATRIX: {   /* This should be replaced with a more general mechanism */
269:     KSPChebyshevSetNewMatrix(ksp);
270:   } break;
271:   default: break;
272:   }

274:   /* scale the matrix if requested */
275:   if (ksp->dscale) {
276:     Mat         mat,pmat;
277:     PetscScalar *xx;
278:     PetscInt    i,n;
279:     PetscBool   zeroflag = PETSC_FALSE;
280:     if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
281:     PCGetOperators(ksp->pc,&mat,&pmat);
282:     if (!ksp->diagonal) { /* allocate vector to hold diagonal */
283:       MatGetVecs(pmat,&ksp->diagonal,0);
284:     }
285:     MatGetDiagonal(pmat,ksp->diagonal);
286:     VecGetLocalSize(ksp->diagonal,&n);
287:     VecGetArray(ksp->diagonal,&xx);
288:     for (i=0; i<n; i++) {
289:       if (xx[i] != 0.0) xx[i] = 1.0/PetscSqrtReal(PetscAbsScalar(xx[i]));
290:       else {
291:         xx[i]    = 1.0;
292:         zeroflag = PETSC_TRUE;
293:       }
294:     }
295:     VecRestoreArray(ksp->diagonal,&xx);
296:     if (zeroflag) {
297:       PetscInfo(ksp,"Zero detected in diagonal of matrix, using 1 at those locations\n");
298:     }
299:     MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
300:     if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
301:     ksp->dscalefix2 = PETSC_FALSE;
302:   }
303:   PetscLogEventEnd(KSP_SetUp,ksp,ksp->vec_rhs,ksp->vec_sol,0);
304:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
305:   PCSetUp(ksp->pc);
306:   if (ksp->nullsp) {
307:     PetscBool test = PETSC_FALSE;
308:     PetscOptionsGetBool(((PetscObject)ksp)->prefix,"-ksp_test_null_space",&test,NULL);
309:     if (test) {
310:       Mat mat;
311:       PCGetOperators(ksp->pc,&mat,NULL);
312:       MatNullSpaceTest(ksp->nullsp,mat,NULL);
313:     }
314:   }
315:   ksp->setupstage = KSP_SETUP_NEWRHS;
316:   return(0);
317: }

319: #include <petscdraw.h>
322: /*@
323:    KSPSolve - Solves linear system.

325:    Collective on KSP

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

332:    Options Database Keys:
333: +  -ksp_compute_eigenvalues - compute preconditioned operators eigenvalues
334: .  -ksp_plot_eigenvalues - plot the computed eigenvalues in an X-window
335: .  -ksp_compute_eigenvalues_explicitly - compute the eigenvalues by forming the dense operator and useing LAPACK
336: .  -ksp_plot_eigenvalues_explicitly - plot the explicitly computing eigenvalues
337: .  -ksp_view_mat binary - save matrix to the default binary viewer
338: .  -ksp_view_pmat binary - save matrix used to build preconditioner to the default binary viewer
339: .  -ksp_view_rhs binary - save right hand side vector to the default binary viewer
340: .  -ksp_view_solution binary - save computed solution vector to the default binary viewer
341:            (can be read later with src/ksp/examples/tutorials/ex10.c for testing solvers)
342: .  -ksp_view_mat_explicit - for matrix-free operators, computes the matrix entries and views them
343: .  -ksp_view_preconditioned_operator_explicit - computes the product of the preconditioner and matrix as an explicit matrix and views it
344: .  -ksp_converged_reason - print reason for converged or diverged, also prints number of iterations
345: .  -ksp_final_residual - print 2-norm of true linear system residual at the end of the solution process
346: -  -ksp_view - print the ksp data structure at the end of the system solution

348:    Notes:

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

352:    The operator is specified with KSPSetOperators().

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

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

362:    Understanding Convergence:
363:    The routines KSPMonitorSet(), KSPComputeEigenvalues(), and
364:    KSPComputeEigenvaluesExplicitly() provide information on additional
365:    options to monitor convergence and print eigenvalue information.

367:    Level: beginner

369: .keywords: KSP, solve, linear system

371: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPConvergedDefault(),
372:           KSPSolveTranspose(), KSPGetIterationNumber()
373: @*/
374: PetscErrorCode  KSPSolve(KSP ksp,Vec b,Vec x)
375: {
376:   PetscErrorCode    ierr;
377:   PetscMPIInt       rank;
378:   PetscBool         flag1,flag2,flag3,flg = PETSC_FALSE,inXisinB=PETSC_FALSE,guess_zero;
379:   Mat               mat,premat;


386:   if (x && x == b) {
387:     if (!ksp->guess_zero) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_INCOMP,"Cannot use x == b with nonzero initial guess");
388:     VecDuplicate(b,&x);
389:     inXisinB = PETSC_TRUE;
390:   }
391:   if (b) {
392:     PetscObjectReference((PetscObject)b);
393:     VecDestroy(&ksp->vec_rhs);
394:     ksp->vec_rhs = b;
395:   }
396:   if (x) {
397:     PetscObjectReference((PetscObject)x);
398:     VecDestroy(&ksp->vec_sol);
399:     ksp->vec_sol = x;
400:   }
401:   KSPViewFromOptions(ksp,NULL,"-ksp_view_pre");

403:   if (ksp->presolve) {
404:     (*ksp->presolve)(ksp,ksp->vec_rhs,ksp->vec_sol,ksp->prectx);
405:   }
406:   PetscLogEventBegin(KSP_Solve,ksp,ksp->vec_rhs,ksp->vec_sol,0);

408:   /* reset the residual history list if requested */
409:   if (ksp->res_hist_reset) ksp->res_hist_len = 0;
410:   ksp->transpose_solve = PETSC_FALSE;

412:   if (ksp->guess) {
413:     KSPFischerGuessFormGuess(ksp->guess,ksp->vec_rhs,ksp->vec_sol);
414:     ksp->guess_zero = PETSC_FALSE;
415:   }
416:   /* KSPSetUp() scales the matrix if needed */
417:   KSPSetUp(ksp);
418:   KSPSetUpOnBlocks(ksp);

420:   /* diagonal scale RHS if called for */
421:   if (ksp->dscale) {
422:     VecPointwiseMult(ksp->vec_rhs,ksp->vec_rhs,ksp->diagonal);
423:     /* second time in, but matrix was scaled back to original */
424:     if (ksp->dscalefix && ksp->dscalefix2) {
425:       Mat mat,pmat;

427:       PCGetOperators(ksp->pc,&mat,&pmat);
428:       MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
429:       if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
430:     }

432:     /*  scale initial guess */
433:     if (!ksp->guess_zero) {
434:       if (!ksp->truediagonal) {
435:         VecDuplicate(ksp->diagonal,&ksp->truediagonal);
436:         VecCopy(ksp->diagonal,ksp->truediagonal);
437:         VecReciprocal(ksp->truediagonal);
438:       }
439:       VecPointwiseMult(ksp->vec_sol,ksp->vec_sol,ksp->truediagonal);
440:     }
441:   }
442:   PCPreSolve(ksp->pc,ksp);

444:   if (ksp->guess_zero) { VecSet(ksp->vec_sol,0.0);}
445:   if (ksp->guess_knoll) {
446:     PCApply(ksp->pc,ksp->vec_rhs,ksp->vec_sol);
447:     KSP_RemoveNullSpace(ksp,ksp->vec_sol);
448:     ksp->guess_zero = PETSC_FALSE;
449:   }

451:   /* can we mark the initial guess as zero for this solve? */
452:   guess_zero = ksp->guess_zero;
453:   if (!ksp->guess_zero) {
454:     PetscReal norm;

456:     VecNormAvailable(ksp->vec_sol,NORM_2,&flg,&norm);
457:     if (flg && !norm) ksp->guess_zero = PETSC_TRUE;
458:   }
459:   (*ksp->ops->solve)(ksp);
460:   ksp->guess_zero = guess_zero;

462:   if (!ksp->reason) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
463:   if (ksp->printreason) {
464:     PetscViewerASCIIAddTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)ksp)),((PetscObject)ksp)->tablevel);
465:     if (ksp->reason > 0) {
466:       PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)ksp)),"Linear solve converged due to %s iterations %D\n",KSPConvergedReasons[ksp->reason],ksp->its);
467:     } else {
468:       PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)ksp)),"Linear solve did not converge due to %s iterations %D\n",KSPConvergedReasons[ksp->reason],ksp->its);
469:     }
470:     PetscViewerASCIISubtractTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)ksp)),((PetscObject)ksp)->tablevel);
471:   }
472:   PCPostSolve(ksp->pc,ksp);

474:   /* diagonal scale solution if called for */
475:   if (ksp->dscale) {
476:     VecPointwiseMult(ksp->vec_sol,ksp->vec_sol,ksp->diagonal);
477:     /* unscale right hand side and matrix */
478:     if (ksp->dscalefix) {
479:       Mat mat,pmat;

481:       VecReciprocal(ksp->diagonal);
482:       VecPointwiseMult(ksp->vec_rhs,ksp->vec_rhs,ksp->diagonal);
483:       PCGetOperators(ksp->pc,&mat,&pmat);
484:       MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
485:       if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
486:       VecReciprocal(ksp->diagonal);
487:       ksp->dscalefix2 = PETSC_TRUE;
488:     }
489:   }
490:   PetscLogEventEnd(KSP_Solve,ksp,ksp->vec_rhs,ksp->vec_sol,0);
491:   if (ksp->postsolve) {
492:     (*ksp->postsolve)(ksp,ksp->vec_rhs,ksp->vec_sol,ksp->postctx);
493:   }

495:   if (ksp->guess) {
496:     KSPFischerGuessUpdate(ksp->guess,ksp->vec_sol);
497:   }

499:   MPI_Comm_rank(PetscObjectComm((PetscObject)ksp),&rank);

501:   PCGetOperators(ksp->pc,&mat,&premat);
502:   MatViewFromOptions(mat,((PetscObject)ksp)->prefix,"-ksp_view_mat");
503:   MatViewFromOptions(premat,((PetscObject)ksp)->prefix,"-ksp_view_pmat");
504:   VecViewFromOptions(ksp->vec_rhs,((PetscObject)ksp)->prefix,"-ksp_view_rhs");

506:   flag1 = PETSC_FALSE;
507:   flag2 = PETSC_FALSE;
508:   flag3 = PETSC_FALSE;
509:   PetscOptionsGetBool(((PetscObject)ksp)->prefix,"-ksp_compute_eigenvalues",&flag1,NULL);
510:   PetscOptionsGetBool(((PetscObject)ksp)->prefix,"-ksp_plot_eigenvalues",&flag2,NULL);
511:   PetscOptionsGetBool(((PetscObject)ksp)->prefix,"-ksp_plot_eigencontours",&flag3,NULL);
512:   if (flag1 || flag2 || flag3) {
513:     PetscInt  nits,n,i,neig;
514:     PetscReal *r,*c;

516:     KSPGetIterationNumber(ksp,&nits);
517:     n    = nits+2;

519:     if (!nits) {
520:       PetscPrintf(PetscObjectComm((PetscObject)ksp),"Zero iterations in solver, cannot approximate any eigenvalues\n");
521:     } else {
522:       PetscMalloc2(n,&r,n,&c);
523:       KSPComputeEigenvalues(ksp,n,r,c,&neig);
524:       if (flag1) {
525:         PetscPrintf(PetscObjectComm((PetscObject)ksp),"Iteratively computed eigenvalues\n");
526:         for (i=0; i<neig; i++) {
527:           if (c[i] >= 0.0) {
528:             PetscPrintf(PetscObjectComm((PetscObject)ksp),"%g + %gi\n",(double)r[i],(double)c[i]);
529:           } else {
530:             PetscPrintf(PetscObjectComm((PetscObject)ksp),"%g - %gi\n",(double)r[i],-(double)c[i]);
531:           }
532:         }
533:       }
534:       if (flag2 && !rank) {
535:         PetscDraw   draw;
536:         PetscDrawSP drawsp;

538:         if (!ksp->eigviewer) {
539:           PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Iteratively Computed Eigenvalues",PETSC_DECIDE,PETSC_DECIDE,400,400,&ksp->eigviewer);
540:         }
541:         PetscViewerDrawGetDraw(ksp->eigviewer,0,&draw);
542:         PetscDrawSPCreate(draw,1,&drawsp);
543:         PetscDrawSPReset(drawsp);
544:         for (i=0; i<neig; i++) {
545:           PetscDrawSPAddPoint(drawsp,r+i,c+i);
546:         }
547:         PetscDrawSPDraw(drawsp,PETSC_TRUE);
548:         PetscDrawSPDestroy(&drawsp);
549:       }
550:       if (flag3 && !rank) {
551:         KSPPlotEigenContours_Private(ksp,neig,r,c);
552:       }
553:       PetscFree2(r,c);
554:     }
555:   }

557:   flag1 = PETSC_FALSE;
558:   PetscOptionsGetBool(((PetscObject)ksp)->prefix,"-ksp_compute_singularvalues",&flag1,NULL);
559:   if (flag1) {
560:     PetscInt nits;

562:     KSPGetIterationNumber(ksp,&nits);
563:     if (!nits) {
564:       PetscPrintf(PetscObjectComm((PetscObject)ksp),"Zero iterations in solver, cannot approximate any singular values\n");
565:     } else {
566:       PetscReal emax,emin;

568:       KSPComputeExtremeSingularValues(ksp,&emax,&emin);
569:       PetscPrintf(PetscObjectComm((PetscObject)ksp),"Iteratively computed extreme singular values: max %g min %g max/min %g\n",(double)emax,(double)emin,(double)(emax/emin));
570:     }
571:   }


574:   flag1 = PETSC_FALSE;
575:   flag2 = PETSC_FALSE;
576:   PetscOptionsGetBool(((PetscObject)ksp)->prefix,"-ksp_compute_eigenvalues_explicitly",&flag1,NULL);
577:   PetscOptionsGetBool(((PetscObject)ksp)->prefix,"-ksp_plot_eigenvalues_explicitly",&flag2,NULL);
578:   if (flag1 || flag2) {
579:     PetscInt  n,i;
580:     PetscReal *r,*c;
581:     VecGetSize(ksp->vec_sol,&n);
582:     PetscMalloc2(n,&r,n,&c);
583:     KSPComputeEigenvaluesExplicitly(ksp,n,r,c);
584:     if (flag1) {
585:       PetscPrintf(PetscObjectComm((PetscObject)ksp),"Explicitly computed eigenvalues\n");
586:       for (i=0; i<n; i++) {
587:         if (c[i] >= 0.0) {
588:           PetscPrintf(PetscObjectComm((PetscObject)ksp),"%g + %gi\n",(double)r[i],(double)c[i]);
589:         } else {
590:           PetscPrintf(PetscObjectComm((PetscObject)ksp),"%g - %gi\n",(double)r[i],-(double)c[i]);
591:         }
592:       }
593:     }
594:     if (flag2 && !rank) {
595:       PetscDraw   draw;
596:       PetscDrawSP drawsp;

598:       if (!ksp->eigviewer) {
599:         PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Explicitly Computed Eigenvalues",0,320,400,400,&ksp->eigviewer);
600:       }
601:       PetscViewerDrawGetDraw(ksp->eigviewer,0,&draw);
602:       PetscDrawSPCreate(draw,1,&drawsp);
603:       PetscDrawSPReset(drawsp);
604:       for (i=0; i<n; i++) {
605:         PetscDrawSPAddPoint(drawsp,r+i,c+i);
606:       }
607:       PetscDrawSPDraw(drawsp,PETSC_TRUE);
608:       PetscDrawSPDestroy(&drawsp);
609:     }
610:     PetscFree2(r,c);
611:   }

613:   PetscOptionsHasName(((PetscObject)ksp)->prefix,"-ksp_view_mat_explicit",&flag2);
614:   if (flag2) {
615:     Mat A,B;
616:     PCGetOperators(ksp->pc,&A,NULL);
617:     MatComputeExplicitOperator(A,&B);
618:     MatViewFromOptions(B,((PetscObject)ksp)->prefix,"-ksp_view_mat_explicit");
619:     MatDestroy(&B);
620:   }
621:   PetscOptionsHasName(((PetscObject)ksp)->prefix,"-ksp_view_preconditioned_operator_explicit",&flag2);
622:   if (flag2) {
623:     Mat B;
624:     KSPComputeExplicitOperator(ksp,&B);
625:     MatViewFromOptions(B,((PetscObject)ksp)->prefix,"-ksp_view_preconditioned_operator_explicit");
626:     MatDestroy(&B);
627:   }
628:   KSPViewFromOptions(ksp,NULL,"-ksp_view");

630:   flg  = PETSC_FALSE;
631:   PetscOptionsGetBool(((PetscObject)ksp)->prefix,"-ksp_final_residual",&flg,NULL);
632:   if (flg) {
633:     Mat       A;
634:     Vec       t;
635:     PetscReal norm;
636:     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");
637:     PCGetOperators(ksp->pc,&A,NULL);
638:     VecDuplicate(ksp->vec_rhs,&t);
639:     KSP_MatMult(ksp,A,ksp->vec_sol,t);
640:     VecAYPX(t, -1.0, ksp->vec_rhs);
641:     VecNorm(t,NORM_2,&norm);
642:     VecDestroy(&t);
643:     PetscPrintf(PetscObjectComm((PetscObject)ksp),"KSP final norm of residual %g\n",(double)norm);
644:   }
645:   VecViewFromOptions(ksp->vec_sol,((PetscObject)ksp)->prefix,"-ksp_view_solution");

647:   if (inXisinB) {
648:     VecCopy(x,b);
649:     VecDestroy(&x);
650:   }
651:   PetscObjectSAWsBlock((PetscObject)ksp);
652:   if (ksp->errorifnotconverged && ksp->reason < 0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged");
653:   return(0);
654: }

658: /*@
659:    KSPSolveTranspose - Solves the transpose of a linear system.

661:    Collective on KSP

663:    Input Parameter:
664: +  ksp - iterative context obtained from KSPCreate()
665: .  b - right hand side vector
666: -  x - solution vector

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

670:    Developer Notes: We need to implement a KSPSolveHermitianTranspose()

672:    Level: developer

674: .keywords: KSP, solve, linear system

676: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPConvergedDefault(),
677:           KSPSolve()
678: @*/

680: PetscErrorCode  KSPSolveTranspose(KSP ksp,Vec b,Vec x)
681: {
683:   PetscBool      inXisinB=PETSC_FALSE;

689:   if (x == b) {
690:     VecDuplicate(b,&x);
691:     inXisinB = PETSC_TRUE;
692:   }
693:   PetscObjectReference((PetscObject)b);
694:   PetscObjectReference((PetscObject)x);
695:   VecDestroy(&ksp->vec_rhs);
696:   VecDestroy(&ksp->vec_sol);

698:   ksp->vec_rhs         = b;
699:   ksp->vec_sol         = x;
700:   ksp->transpose_solve = PETSC_TRUE;

702:   KSPSetUp(ksp);
703:   if (ksp->guess_zero) { VecSet(ksp->vec_sol,0.0);}
704:   (*ksp->ops->solve)(ksp);
705:   if (!ksp->reason) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
706:   if (ksp->printreason) {
707:     if (ksp->reason > 0) {
708:       PetscPrintf(PetscObjectComm((PetscObject)ksp),"Linear solve converged due to %s iterations %D\n",KSPConvergedReasons[ksp->reason],ksp->its);
709:     } else {
710:       PetscPrintf(PetscObjectComm((PetscObject)ksp),"Linear solve did not converge due to %s iterations %D\n",KSPConvergedReasons[ksp->reason],ksp->its);
711:     }
712:   }
713:   if (inXisinB) {
714:     VecCopy(x,b);
715:     VecDestroy(&x);
716:   }
717:   if (ksp->errorifnotconverged && ksp->reason < 0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged");
718:   return(0);
719: }

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

726:    Collective on KSP

728:    Input Parameter:
729: .  ksp - iterative context obtained from KSPCreate()

731:    Level: beginner

733: .keywords: KSP, destroy

735: .seealso: KSPCreate(), KSPSetUp(), KSPSolve()
736: @*/
737: PetscErrorCode  KSPReset(KSP ksp)
738: {

743:   if (!ksp) return(0);
744:   if (ksp->ops->reset) {
745:     (*ksp->ops->reset)(ksp);
746:   }
747:   if (ksp->pc) {PCReset(ksp->pc);}
748:   KSPFischerGuessDestroy(&ksp->guess);
749:   VecDestroyVecs(ksp->nwork,&ksp->work);
750:   VecDestroy(&ksp->vec_rhs);
751:   VecDestroy(&ksp->vec_sol);
752:   VecDestroy(&ksp->diagonal);
753:   VecDestroy(&ksp->truediagonal);
754:   MatNullSpaceDestroy(&ksp->nullsp);

756:   ksp->setupstage = KSP_SETUP_NEW;
757:   return(0);
758: }

762: /*@
763:    KSPDestroy - Destroys KSP context.

765:    Collective on KSP

767:    Input Parameter:
768: .  ksp - iterative context obtained from KSPCreate()

770:    Level: beginner

772: .keywords: KSP, destroy

774: .seealso: KSPCreate(), KSPSetUp(), KSPSolve()
775: @*/
776: PetscErrorCode  KSPDestroy(KSP *ksp)
777: {
779:   PC             pc;

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

786:   PetscObjectSAWsViewOff((PetscObject)*ksp);
787:   /*
788:    Avoid a cascading call to PCReset(ksp->pc) from the following call:
789:    PCReset() shouldn't be called from KSPDestroy() as it is unprotected by pc's
790:    refcount (and may be shared, e.g., by other ksps).
791:    */
792:   pc         = (*ksp)->pc;
793:   (*ksp)->pc = NULL;
794:   KSPReset((*ksp));
795:   (*ksp)->pc = pc;
796:     if ((*ksp)->ops->destroy) {(*(*ksp)->ops->destroy)(*ksp);}

798:   DMDestroy(&(*ksp)->dm);
799:   PCDestroy(&(*ksp)->pc);
800:   PetscFree((*ksp)->res_hist_alloc);
801:   if ((*ksp)->convergeddestroy) {
802:     (*(*ksp)->convergeddestroy)((*ksp)->cnvP);
803:   }
804:   KSPMonitorCancel((*ksp));
805:   PetscViewerDestroy(&(*ksp)->eigviewer);
806:   PetscHeaderDestroy(ksp);
807:   return(0);
808: }

812: /*@
813:     KSPSetPCSide - Sets the preconditioning side.

815:     Logically Collective on KSP

817:     Input Parameter:
818: .   ksp - iterative context obtained from KSPCreate()

820:     Output Parameter:
821: .   side - the preconditioning side, where side is one of
822: .vb
823:       PC_LEFT - left preconditioning (default)
824:       PC_RIGHT - right preconditioning
825:       PC_SYMMETRIC - symmetric preconditioning
826: .ve

828:     Options Database Keys:
829: .   -ksp_pc_side <right,left,symmetric>

831:     Notes:
832:     Left preconditioning is used by default for most Krylov methods except KSPFGMRES which only supports right preconditioning.
833:     Symmetric preconditioning is currently available only for the KSPQCG method. Note, however, that
834:     symmetric preconditioning can be emulated by using either right or left
835:     preconditioning and a pre or post processing step.

837:     Level: intermediate

839: .keywords: KSP, set, right, left, symmetric, side, preconditioner, flag

841: .seealso: KSPGetPCSide()
842: @*/
843: PetscErrorCode  KSPSetPCSide(KSP ksp,PCSide side)
844: {
848:   ksp->pc_side = ksp->pc_side_set = side;
849:   return(0);
850: }

854: /*@
855:     KSPGetPCSide - Gets the preconditioning side.

857:     Not Collective

859:     Input Parameter:
860: .   ksp - iterative context obtained from KSPCreate()

862:     Output Parameter:
863: .   side - the preconditioning side, where side is one of
864: .vb
865:       PC_LEFT - left preconditioning (default)
866:       PC_RIGHT - right preconditioning
867:       PC_SYMMETRIC - symmetric preconditioning
868: .ve

870:     Level: intermediate

872: .keywords: KSP, get, right, left, symmetric, side, preconditioner, flag

874: .seealso: KSPSetPCSide()
875: @*/
876: PetscErrorCode  KSPGetPCSide(KSP ksp,PCSide *side)
877: {

883:   KSPSetUpNorms_Private(ksp,&ksp->normtype,&ksp->pc_side);
884:   *side = ksp->pc_side;
885:   return(0);
886: }

890: /*@
891:    KSPGetTolerances - Gets the relative, absolute, divergence, and maximum
892:    iteration tolerances used by the default KSP convergence tests.

894:    Not Collective

896:    Input Parameter:
897: .  ksp - the Krylov subspace context

899:    Output Parameters:
900: +  rtol - the relative convergence tolerance
901: .  abstol - the absolute convergence tolerance
902: .  dtol - the divergence tolerance
903: -  maxits - maximum number of iterations

905:    Notes:
906:    The user can specify NULL for any parameter that is not needed.

908:    Level: intermediate

910: .keywords: KSP, get, tolerance, absolute, relative, divergence, convergence,
911:            maximum, iterations

913: .seealso: KSPSetTolerances()
914: @*/
915: PetscErrorCode  KSPGetTolerances(KSP ksp,PetscReal *rtol,PetscReal *abstol,PetscReal *dtol,PetscInt *maxits)
916: {
919:   if (abstol) *abstol = ksp->abstol;
920:   if (rtol) *rtol = ksp->rtol;
921:   if (dtol) *dtol = ksp->divtol;
922:   if (maxits) *maxits = ksp->max_it;
923:   return(0);
924: }

928: /*@
929:    KSPSetTolerances - Sets the relative, absolute, divergence, and maximum
930:    iteration tolerances used by the default KSP convergence testers.

932:    Logically Collective on KSP

934:    Input Parameters:
935: +  ksp - the Krylov subspace context
936: .  rtol - the relative convergence tolerance, relative decrease in the (possibly preconditioned) residual norm
937: .  abstol - the absolute convergence tolerance   absolute size of the (possibly preconditioned) residual norm
938: .  dtol - the divergence tolerance,   amount (possibly preconditioned) residual norm can increase before KSPConvergedDefault() concludes that the method is diverging
939: -  maxits - maximum number of iterations to use

941:    Options Database Keys:
942: +  -ksp_atol <abstol> - Sets abstol
943: .  -ksp_rtol <rtol> - Sets rtol
944: .  -ksp_divtol <dtol> - Sets dtol
945: -  -ksp_max_it <maxits> - Sets maxits

947:    Notes:
948:    Use PETSC_DEFAULT to retain the default value of any of the tolerances.

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

953:    Level: intermediate

955: .keywords: KSP, set, tolerance, absolute, relative, divergence,
956:            convergence, maximum, iterations

958: .seealso: KSPGetTolerances(), KSPConvergedDefault(), KSPSetConvergenceTest()
959: @*/
960: PetscErrorCode  KSPSetTolerances(KSP ksp,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt maxits)
961: {

969:   if (rtol != PETSC_DEFAULT) {
970:     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);
971:     ksp->rtol = rtol;
972:   }
973:   if (abstol != PETSC_DEFAULT) {
974:     if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
975:     ksp->abstol = abstol;
976:   }
977:   if (dtol != PETSC_DEFAULT) {
978:     if (dtol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Divergence tolerance %g must be larger than 1.0",(double)dtol);
979:     ksp->divtol = dtol;
980:   }
981:   if (maxits != PETSC_DEFAULT) {
982:     if (maxits < 0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxits);
983:     ksp->max_it = maxits;
984:   }
985:   return(0);
986: }

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

995:    Logically Collective on KSP

997:    Input Parameters:
998: +  ksp - iterative context obtained from KSPCreate()
999: -  flg - PETSC_TRUE indicates the guess is non-zero, PETSC_FALSE indicates the guess is zero

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

1004:    Level: beginner

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

1009: .keywords: KSP, set, initial guess, nonzero

1011: .seealso: KSPGetInitialGuessNonzero(), KSPSetInitialGuessKnoll(), KSPGetInitialGuessKnoll()
1012: @*/
1013: PetscErrorCode  KSPSetInitialGuessNonzero(KSP ksp,PetscBool flg)
1014: {
1018:   ksp->guess_zero = (PetscBool) !(int)flg;
1019:   return(0);
1020: }

1024: /*@
1025:    KSPGetInitialGuessNonzero - Determines whether the KSP solver is using
1026:    a zero initial guess.

1028:    Not Collective

1030:    Input Parameter:
1031: .  ksp - iterative context obtained from KSPCreate()

1033:    Output Parameter:
1034: .  flag - PETSC_TRUE if guess is nonzero, else PETSC_FALSE

1036:    Level: intermediate

1038: .keywords: KSP, set, initial guess, nonzero

1040: .seealso: KSPSetInitialGuessNonzero(), KSPSetInitialGuessKnoll(), KSPGetInitialGuessKnoll()
1041: @*/
1042: PetscErrorCode  KSPGetInitialGuessNonzero(KSP ksp,PetscBool  *flag)
1043: {
1047:   if (ksp->guess_zero) *flag = PETSC_FALSE;
1048:   else *flag = PETSC_TRUE;
1049:   return(0);
1050: }

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

1057:    Logically Collective on KSP

1059:    Input Parameters:
1060: +  ksp - iterative context obtained from KSPCreate()
1061: -  flg - PETSC_TRUE indicates you want the error generated

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

1066:    Level: intermediate

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

1072: .keywords: KSP, set, initial guess, nonzero

1074: .seealso: KSPGetInitialGuessNonzero(), KSPSetInitialGuessKnoll(), KSPGetInitialGuessKnoll(), KSPGetErrorIfNotConverged()
1075: @*/
1076: PetscErrorCode  KSPSetErrorIfNotConverged(KSP ksp,PetscBool flg)
1077: {
1081:   ksp->errorifnotconverged = flg;
1082:   return(0);
1083: }

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

1090:    Not Collective

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

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

1098:    Level: intermediate

1100: .keywords: KSP, set, initial guess, nonzero

1102: .seealso: KSPSetInitialGuessNonzero(), KSPSetInitialGuessKnoll(), KSPGetInitialGuessKnoll(), KSPSetErrorIfNotConverged()
1103: @*/
1104: PetscErrorCode  KSPGetErrorIfNotConverged(KSP ksp,PetscBool  *flag)
1105: {
1109:   *flag = ksp->errorifnotconverged;
1110:   return(0);
1111: }

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

1118:    Logically Collective on KSP

1120:    Input Parameters:
1121: +  ksp - iterative context obtained from KSPCreate()
1122: -  flg - PETSC_TRUE or PETSC_FALSE

1124:    Level: advanced


1127: .keywords: KSP, set, initial guess, nonzero

1129: .seealso: KSPGetInitialGuessKnoll(), KSPSetInitialGuessNonzero(), KSPGetInitialGuessNonzero()
1130: @*/
1131: PetscErrorCode  KSPSetInitialGuessKnoll(KSP ksp,PetscBool flg)
1132: {
1136:   ksp->guess_knoll = flg;
1137:   return(0);
1138: }

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

1146:    Not Collective

1148:    Input Parameter:
1149: .  ksp - iterative context obtained from KSPCreate()

1151:    Output Parameter:
1152: .  flag - PETSC_TRUE if using Knoll trick, else PETSC_FALSE

1154:    Level: advanced

1156: .keywords: KSP, set, initial guess, nonzero

1158: .seealso: KSPSetInitialGuessKnoll(), KSPSetInitialGuessNonzero(), KSPGetInitialGuessNonzero()
1159: @*/
1160: PetscErrorCode  KSPGetInitialGuessKnoll(KSP ksp,PetscBool  *flag)
1161: {
1165:   *flag = ksp->guess_knoll;
1166:   return(0);
1167: }

1171: /*@
1172:    KSPGetComputeSingularValues - Gets the flag indicating whether the extreme singular
1173:    values will be calculated via a Lanczos or Arnoldi process as the linear
1174:    system is solved.

1176:    Not Collective

1178:    Input Parameter:
1179: .  ksp - iterative context obtained from KSPCreate()

1181:    Output Parameter:
1182: .  flg - PETSC_TRUE or PETSC_FALSE

1184:    Options Database Key:
1185: .  -ksp_monitor_singular_value - Activates KSPSetComputeSingularValues()

1187:    Notes:
1188:    Currently this option is not valid for all iterative methods.

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

1194:    Level: advanced

1196: .keywords: KSP, set, compute, singular values

1198: .seealso: KSPComputeExtremeSingularValues(), KSPMonitorSingularValue()
1199: @*/
1200: PetscErrorCode  KSPGetComputeSingularValues(KSP ksp,PetscBool  *flg)
1201: {
1205:   *flg = ksp->calc_sings;
1206:   return(0);
1207: }

1211: /*@
1212:    KSPSetComputeSingularValues - Sets a flag so that the extreme singular
1213:    values will be calculated via a Lanczos or Arnoldi process as the linear
1214:    system is solved.

1216:    Logically Collective on KSP

1218:    Input Parameters:
1219: +  ksp - iterative context obtained from KSPCreate()
1220: -  flg - PETSC_TRUE or PETSC_FALSE

1222:    Options Database Key:
1223: .  -ksp_monitor_singular_value - Activates KSPSetComputeSingularValues()

1225:    Notes:
1226:    Currently this option is not valid for all iterative methods.

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

1232:    Level: advanced

1234: .keywords: KSP, set, compute, singular values

1236: .seealso: KSPComputeExtremeSingularValues(), KSPMonitorSingularValue()
1237: @*/
1238: PetscErrorCode  KSPSetComputeSingularValues(KSP ksp,PetscBool flg)
1239: {
1243:   ksp->calc_sings = flg;
1244:   return(0);
1245: }

1249: /*@
1250:    KSPGetComputeEigenvalues - Gets the flag indicating that the extreme eigenvalues
1251:    values will be calculated via a Lanczos or Arnoldi process as the linear
1252:    system is solved.

1254:    Not Collective

1256:    Input Parameter:
1257: .  ksp - iterative context obtained from KSPCreate()

1259:    Output Parameter:
1260: .  flg - PETSC_TRUE or PETSC_FALSE

1262:    Notes:
1263:    Currently this option is not valid for all iterative methods.

1265:    Level: advanced

1267: .keywords: KSP, set, compute, eigenvalues

1269: .seealso: KSPComputeEigenvalues(), KSPComputeEigenvaluesExplicitly()
1270: @*/
1271: PetscErrorCode  KSPGetComputeEigenvalues(KSP ksp,PetscBool  *flg)
1272: {
1276:   *flg = ksp->calc_sings;
1277:   return(0);
1278: }

1282: /*@
1283:    KSPSetComputeEigenvalues - Sets a flag so that the extreme eigenvalues
1284:    values will be calculated via a Lanczos or Arnoldi process as the linear
1285:    system is solved.

1287:    Logically Collective on KSP

1289:    Input Parameters:
1290: +  ksp - iterative context obtained from KSPCreate()
1291: -  flg - PETSC_TRUE or PETSC_FALSE

1293:    Notes:
1294:    Currently this option is not valid for all iterative methods.

1296:    Level: advanced

1298: .keywords: KSP, set, compute, eigenvalues

1300: .seealso: KSPComputeEigenvalues(), KSPComputeEigenvaluesExplicitly()
1301: @*/
1302: PetscErrorCode  KSPSetComputeEigenvalues(KSP ksp,PetscBool flg)
1303: {
1307:   ksp->calc_sings = flg;
1308:   return(0);
1309: }

1313: /*@
1314:    KSPGetRhs - Gets the right-hand-side vector for the linear system to
1315:    be solved.

1317:    Not Collective

1319:    Input Parameter:
1320: .  ksp - iterative context obtained from KSPCreate()

1322:    Output Parameter:
1323: .  r - right-hand-side vector

1325:    Level: developer

1327: .keywords: KSP, get, right-hand-side, rhs

1329: .seealso: KSPGetSolution(), KSPSolve()
1330: @*/
1331: PetscErrorCode  KSPGetRhs(KSP ksp,Vec *r)
1332: {
1336:   *r = ksp->vec_rhs;
1337:   return(0);
1338: }

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

1347:    Not Collective

1349:    Input Parameters:
1350: .  ksp - iterative context obtained from KSPCreate()

1352:    Output Parameters:
1353: .  v - solution vector

1355:    Level: developer

1357: .keywords: KSP, get, solution

1359: .seealso: KSPGetRhs(),  KSPBuildSolution(), KSPSolve()
1360: @*/
1361: PetscErrorCode  KSPGetSolution(KSP ksp,Vec *v)
1362: {
1366:   *v = ksp->vec_sol;
1367:   return(0);
1368: }

1372: /*@
1373:    KSPSetPC - Sets the preconditioner to be used to calculate the
1374:    application of the preconditioner on a vector.

1376:    Collective on KSP

1378:    Input Parameters:
1379: +  ksp - iterative context obtained from KSPCreate()
1380: -  pc   - the preconditioner object

1382:    Notes:
1383:    Use KSPGetPC() to retrieve the preconditioner context (for example,
1384:    to free it at the end of the computations).

1386:    Level: developer

1388: .keywords: KSP, set, precondition, Binv

1390: .seealso: KSPGetPC()
1391: @*/
1392: PetscErrorCode  KSPSetPC(KSP ksp,PC pc)
1393: {

1400:   PetscObjectReference((PetscObject)pc);
1401:   PCDestroy(&ksp->pc);
1402:   ksp->pc = pc;
1403:   PetscLogObjectParent((PetscObject)ksp,(PetscObject)ksp->pc);
1404:   return(0);
1405: }

1409: /*@
1410:    KSPGetPC - Returns a pointer to the preconditioner context
1411:    set with KSPSetPC().

1413:    Not Collective

1415:    Input Parameters:
1416: .  ksp - iterative context obtained from KSPCreate()

1418:    Output Parameter:
1419: .  pc - preconditioner context

1421:    Level: developer

1423: .keywords: KSP, get, preconditioner, Binv

1425: .seealso: KSPSetPC()
1426: @*/
1427: PetscErrorCode  KSPGetPC(KSP ksp,PC *pc)
1428: {

1434:   if (!ksp->pc) {
1435:     PCCreate(PetscObjectComm((PetscObject)ksp),&ksp->pc);
1436:     PetscObjectIncrementTabLevel((PetscObject)ksp->pc,(PetscObject)ksp,0);
1437:     PetscLogObjectParent((PetscObject)ksp,(PetscObject)ksp->pc);
1438:   }
1439:   *pc = ksp->pc;
1440:   return(0);
1441: }

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

1448:    Collective on KSP

1450:    Input Parameters:
1451: +  ksp - iterative context obtained from KSPCreate()
1452: .  it - iteration number
1453: -  rnorm - relative norm of the residual

1455:    Notes:
1456:    This routine is called by the KSP implementations.
1457:    It does not typically need to be called by the user.

1459:    Level: developer

1461: .seealso: KSPMonitorSet()
1462: @*/
1463: PetscErrorCode KSPMonitor(KSP ksp,PetscInt it,PetscReal rnorm)
1464: {
1465:   PetscInt       i, n = ksp->numbermonitors;

1469:   for (i=0; i<n; i++) {
1470:     (*ksp->monitor[i])(ksp,it,rnorm,ksp->monitorcontext[i]);
1471:   }
1472:   return(0);
1473: }

1477: /*@C
1478:    KSPMonitorSet - Sets an ADDITIONAL function to be called at every iteration to monitor
1479:    the residual/error etc.

1481:    Logically Collective on KSP

1483:    Input Parameters:
1484: +  ksp - iterative context obtained from KSPCreate()
1485: .  monitor - pointer to function (if this is NULL, it turns off monitoring
1486: .  mctx    - [optional] context for private data for the
1487:              monitor routine (use NULL if no context is desired)
1488: -  monitordestroy - [optional] routine that frees monitor context
1489:           (may be NULL)

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

1494: +  ksp - iterative context obtained from KSPCreate()
1495: .  it - iteration number
1496: .  rnorm - (estimated) 2-norm of (preconditioned) residual
1497: -  mctx  - optional monitoring context, as set by KSPMonitorSet()

1499:    Options Database Keys:
1500: +    -ksp_monitor        - sets KSPMonitorDefault()
1501: .    -ksp_monitor_true_residual    - sets KSPMonitorTrueResidualNorm()
1502: .    -ksp_monitor_max    - sets KSPMonitorTrueResidualMaxNorm()
1503: .    -ksp_monitor_lg_residualnorm    - sets line graph monitor,
1504:                            uses KSPMonitorLGResidualNormCreate()
1505: .    -ksp_monitor_lg_true_residualnorm   - sets line graph monitor,
1506:                            uses KSPMonitorLGResidualNormCreate()
1507: .    -ksp_monitor_singular_value    - sets KSPMonitorSingularValue()
1508: -    -ksp_monitor_cancel - cancels all monitors that have
1509:                           been hardwired into a code by
1510:                           calls to KSPMonitorSet(), but
1511:                           does not cancel those set via
1512:                           the options database.

1514:    Notes:
1515:    The default is to do nothing.  To print the residual, or preconditioned
1516:    residual if KSPSetNormType(ksp,KSP_NORM_PRECONDITIONED) was called, use
1517:    KSPMonitorDefault() as the monitoring routine, with a null monitoring
1518:    context.

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

1524:    Fortran notes: Only a single monitor function can be set for each KSP object

1526:    Level: beginner

1528: .keywords: KSP, set, monitor

1530: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(), KSPMonitorCancel()
1531: @*/
1532: PetscErrorCode  KSPMonitorSet(KSP ksp,PetscErrorCode (*monitor)(KSP,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
1533: {
1534:   PetscInt       i;

1539:   if (ksp->numbermonitors >= MAXKSPMONITORS) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Too many KSP monitors set");
1540:   for (i=0; i<ksp->numbermonitors;i++) {
1541:     if (monitor == ksp->monitor[i] && monitordestroy == ksp->monitordestroy[i] && mctx == ksp->monitorcontext[i]) {
1542:       if (monitordestroy) {
1543:         (*monitordestroy)(&mctx);
1544:       }
1545:       return(0);
1546:     }
1547:   }
1548:   ksp->monitor[ksp->numbermonitors]          = monitor;
1549:   ksp->monitordestroy[ksp->numbermonitors]   = monitordestroy;
1550:   ksp->monitorcontext[ksp->numbermonitors++] = (void*)mctx;
1551:   return(0);
1552: }

1556: /*@
1557:    KSPMonitorCancel - Clears all monitors for a KSP object.

1559:    Logically Collective on KSP

1561:    Input Parameters:
1562: .  ksp - iterative context obtained from KSPCreate()

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

1569:    Level: intermediate

1571: .keywords: KSP, set, monitor

1573: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(), KSPMonitorSet()
1574: @*/
1575: PetscErrorCode  KSPMonitorCancel(KSP ksp)
1576: {
1578:   PetscInt       i;

1582:   for (i=0; i<ksp->numbermonitors; i++) {
1583:     if (ksp->monitordestroy[i]) {
1584:       (*ksp->monitordestroy[i])(&ksp->monitorcontext[i]);
1585:     }
1586:   }
1587:   ksp->numbermonitors = 0;
1588:   return(0);
1589: }

1593: /*@C
1594:    KSPGetMonitorContext - Gets the monitoring context, as set by
1595:    KSPMonitorSet() for the FIRST monitor only.

1597:    Not Collective

1599:    Input Parameter:
1600: .  ksp - iterative context obtained from KSPCreate()

1602:    Output Parameter:
1603: .  ctx - monitoring context

1605:    Level: intermediate

1607: .keywords: KSP, get, monitor, context

1609: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate()
1610: @*/
1611: PetscErrorCode  KSPGetMonitorContext(KSP ksp,void **ctx)
1612: {
1615:   *ctx =      (ksp->monitorcontext[0]);
1616:   return(0);
1617: }

1621: /*@
1622:    KSPSetResidualHistory - Sets the array used to hold the residual history.
1623:    If set, this array will contain the residual norms computed at each
1624:    iteration of the solver.

1626:    Not Collective

1628:    Input Parameters:
1629: +  ksp - iterative context obtained from KSPCreate()
1630: .  a   - array to hold history
1631: .  na  - size of a
1632: -  reset - PETSC_TRUE indicates the history counter is reset to zero
1633:            for each new linear solve

1635:    Level: advanced

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

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

1643: .keywords: KSP, set, residual, history, norm

1645: .seealso: KSPGetResidualHistory()

1647: @*/
1648: PetscErrorCode  KSPSetResidualHistory(KSP ksp,PetscReal a[],PetscInt na,PetscBool reset)
1649: {


1655:   PetscFree(ksp->res_hist_alloc);
1656:   if (na != PETSC_DECIDE && na != PETSC_DEFAULT && a) {
1657:     ksp->res_hist     = a;
1658:     ksp->res_hist_max = na;
1659:   } else {
1660:     if (na != PETSC_DECIDE && na != PETSC_DEFAULT) ksp->res_hist_max = na;
1661:     else                                           ksp->res_hist_max = 10000; /* like default ksp->max_it */
1662:     PetscCalloc1(ksp->res_hist_max,&ksp->res_hist_alloc);

1664:     ksp->res_hist = ksp->res_hist_alloc;
1665:   }
1666:   ksp->res_hist_len   = 0;
1667:   ksp->res_hist_reset = reset;
1668:   return(0);
1669: }

1673: /*@C
1674:    KSPGetResidualHistory - Gets the array used to hold the residual history
1675:    and the number of residuals it contains.

1677:    Not Collective

1679:    Input Parameter:
1680: .  ksp - iterative context obtained from KSPCreate()

1682:    Output Parameters:
1683: +  a   - pointer to array to hold history (or NULL)
1684: -  na  - number of used entries in a (or NULL)

1686:    Level: advanced

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

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

1697: .keywords: KSP, get, residual, history, norm

1699: .seealso: KSPGetResidualHistory()

1701: @*/
1702: PetscErrorCode  KSPGetResidualHistory(KSP ksp,PetscReal *a[],PetscInt *na)
1703: {
1706:   if (a) *a = ksp->res_hist;
1707:   if (na) *na = ksp->res_hist_len;
1708:   return(0);
1709: }

1713: /*@C
1714:    KSPSetConvergenceTest - Sets the function to be used to determine
1715:    convergence.

1717:    Logically Collective on KSP

1719:    Input Parameters:
1720: +  ksp - iterative context obtained from KSPCreate()
1721: .  converge - pointer to int function
1722: .  cctx    - context for private data for the convergence routine (may be null)
1723: -  destroy - a routine for destroying the context (may be null)

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

1728: +  ksp - iterative context obtained from KSPCreate()
1729: .  it - iteration number
1730: .  rnorm - (estimated) 2-norm of (preconditioned) residual
1731: .  reason - the reason why it has converged or diverged
1732: -  cctx  - optional convergence context, as set by KSPSetConvergenceTest()


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

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

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

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

1749:    Level: advanced

1751: .keywords: KSP, set, convergence, test, context

1753: .seealso: KSPConvergedDefault(), KSPGetConvergenceContext(), KSPSetTolerances()
1754: @*/
1755: PetscErrorCode  KSPSetConvergenceTest(KSP ksp,PetscErrorCode (*converge)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
1756: {

1761:   if (ksp->convergeddestroy) {
1762:     (*ksp->convergeddestroy)(ksp->cnvP);
1763:   }
1764:   ksp->converged        = converge;
1765:   ksp->convergeddestroy = destroy;
1766:   ksp->cnvP             = (void*)cctx;
1767:   return(0);
1768: }

1772: /*@C
1773:    KSPGetConvergenceContext - Gets the convergence context set with
1774:    KSPSetConvergenceTest().

1776:    Not Collective

1778:    Input Parameter:
1779: .  ksp - iterative context obtained from KSPCreate()

1781:    Output Parameter:
1782: .  ctx - monitoring context

1784:    Level: advanced

1786: .keywords: KSP, get, convergence, test, context

1788: .seealso: KSPConvergedDefault(), KSPSetConvergenceTest()
1789: @*/
1790: PetscErrorCode  KSPGetConvergenceContext(KSP ksp,void **ctx)
1791: {
1794:   *ctx = ksp->cnvP;
1795:   return(0);
1796: }

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

1804:    Collective on KSP

1806:    Input Parameter:
1807: .  ctx - iterative context obtained from KSPCreate()

1809:    Output Parameter:
1810:    Provide exactly one of
1811: +  v - location to stash solution.
1812: -  V - the solution is returned in this location. This vector is created
1813:        internally. This vector should NOT be destroyed by the user with
1814:        VecDestroy().

1816:    Notes:
1817:    This routine can be used in one of two ways
1818: .vb
1819:       KSPBuildSolution(ksp,NULL,&V);
1820:    or
1821:       KSPBuildSolution(ksp,v,NULL); or KSPBuildSolution(ksp,v,&v);
1822: .ve
1823:    In the first case an internal vector is allocated to store the solution
1824:    (the user cannot destroy this vector). In the second case the solution
1825:    is generated in the vector that the user provides. Note that for certain
1826:    methods, such as KSPCG, the second case requires a copy of the solution,
1827:    while in the first case the call is essentially free since it simply
1828:    returns the vector where the solution already is stored. For some methods
1829:    like GMRES this is a reasonably expensive operation and should only be
1830:    used in truly needed.

1832:    Level: advanced

1834: .keywords: KSP, build, solution

1836: .seealso: KSPGetSolution(), KSPBuildResidual()
1837: @*/
1838: PetscErrorCode  KSPBuildSolution(KSP ksp,Vec v,Vec *V)
1839: {

1844:   if (!V && !v) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONG,"Must provide either v or V");
1845:   if (!V) V = &v;
1846:   (*ksp->ops->buildsolution)(ksp,v,V);
1847:   return(0);
1848: }

1852: /*@C
1853:    KSPBuildResidual - Builds the residual in a vector provided.

1855:    Collective on KSP

1857:    Input Parameter:
1858: .  ksp - iterative context obtained from KSPCreate()

1860:    Output Parameters:
1861: +  v - optional location to stash residual.  If v is not provided,
1862:        then a location is generated.
1863: .  t - work vector.  If not provided then one is generated.
1864: -  V - the residual

1866:    Notes:
1867:    Regardless of whether or not v is provided, the residual is
1868:    returned in V.

1870:    Level: advanced

1872: .keywords: KSP, build, residual

1874: .seealso: KSPBuildSolution()
1875: @*/
1876: PetscErrorCode  KSPBuildResidual(KSP ksp,Vec t,Vec v,Vec *V)
1877: {
1879:   PetscBool      flag = PETSC_FALSE;
1880:   Vec            w    = v,tt = t;

1884:   if (!w) {
1885:     VecDuplicate(ksp->vec_rhs,&w);
1886:     PetscLogObjectParent((PetscObject)ksp,(PetscObject)w);
1887:   }
1888:   if (!tt) {
1889:     VecDuplicate(ksp->vec_sol,&tt); flag = PETSC_TRUE;
1890:     PetscLogObjectParent((PetscObject)ksp,(PetscObject)tt);
1891:   }
1892:   (*ksp->ops->buildresidual)(ksp,tt,w,V);
1893:   if (flag) {VecDestroy(&tt);}
1894:   return(0);
1895: }

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

1903:    Logically Collective on KSP

1905:    Input Parameter:
1906: +  ksp - the KSP context
1907: -  scale - PETSC_TRUE or PETSC_FALSE

1909:    Options Database Key:
1910: +   -ksp_diagonal_scale -
1911: -   -ksp_diagonal_scale_fix - scale the matrix back AFTER the solve


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

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

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

1924:     If you use this with the PCType Eisenstat preconditioner than you can
1925:     use the PCEisenstatNoDiagonalScaling() option, or -pc_eisenstat_no_diagonal_scaling
1926:     to save some unneeded, redundant flops.

1928:    Level: intermediate

1930: .keywords: KSP, set, options, prefix, database

1932: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScaleFix()
1933: @*/
1934: PetscErrorCode  KSPSetDiagonalScale(KSP ksp,PetscBool scale)
1935: {
1939:   ksp->dscale = scale;
1940:   return(0);
1941: }

1945: /*@
1946:    KSPGetDiagonalScale - Checks if KSP solver scales the matrix and
1947:                           right hand side

1949:    Not Collective

1951:    Input Parameter:
1952: .  ksp - the KSP context

1954:    Output Parameter:
1955: .  scale - PETSC_TRUE or PETSC_FALSE

1957:    Notes:
1958:     BE CAREFUL with this routine: it actually scales the matrix and right
1959:     hand side that define the system. After the system is solved the matrix
1960:     and right hand side remain scaled  unless you use KSPSetDiagonalScaleFix()

1962:    Level: intermediate

1964: .keywords: KSP, set, options, prefix, database

1966: .seealso: KSPSetDiagonalScale(), KSPSetDiagonalScaleFix()
1967: @*/
1968: PetscErrorCode  KSPGetDiagonalScale(KSP ksp,PetscBool  *scale)
1969: {
1973:   *scale = ksp->dscale;
1974:   return(0);
1975: }

1979: /*@
1980:    KSPSetDiagonalScaleFix - Tells KSP to diagonally scale the system
1981:      back after solving.

1983:    Logically Collective on KSP

1985:    Input Parameter:
1986: +  ksp - the KSP context
1987: -  fix - PETSC_TRUE to scale back after the system solve, PETSC_FALSE to not
1988:          rescale (default)

1990:    Notes:
1991:      Must be called after KSPSetDiagonalScale()

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

1998:    Level: intermediate

2000: .keywords: KSP, set, options, prefix, database

2002: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScale(), KSPGetDiagonalScaleFix()
2003: @*/
2004: PetscErrorCode  KSPSetDiagonalScaleFix(KSP ksp,PetscBool fix)
2005: {
2009:   ksp->dscalefix = fix;
2010:   return(0);
2011: }

2015: /*@
2016:    KSPGetDiagonalScaleFix - Determines if KSP diagonally scales the system
2017:      back after solving.

2019:    Not Collective

2021:    Input Parameter:
2022: .  ksp - the KSP context

2024:    Output Parameter:
2025: .  fix - PETSC_TRUE to scale back after the system solve, PETSC_FALSE to not
2026:          rescale (default)

2028:    Notes:
2029:      Must be called after KSPSetDiagonalScale()

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

2036:    Level: intermediate

2038: .keywords: KSP, set, options, prefix, database

2040: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScale(), KSPSetDiagonalScaleFix()
2041: @*/
2042: PetscErrorCode  KSPGetDiagonalScaleFix(KSP ksp,PetscBool  *fix)
2043: {
2047:   *fix = ksp->dscalefix;
2048:   return(0);
2049: }

2053: /*@C
2054:    KSPSetComputeOperators - set routine to compute the linear operators

2056:    Logically Collective

2058:    Input Arguments:
2059: +  ksp - the KSP context
2060: .  func - function to compute the operators
2061: -  ctx - optional context

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

2066: +  ksp - the KSP context
2067: .  A - the linear operator
2068: .  B - preconditioning matrix
2069: -  ctx - optional user-provided context

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

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

2076:    Level: beginner

2078: .seealso: KSPSetOperators(), KSPSetComputeRHS(), DMKSPSetComputeOperators(), KSPSetComputeInitialGuess()
2079: @*/
2080: PetscErrorCode KSPSetComputeOperators(KSP ksp,PetscErrorCode (*func)(KSP,Mat,Mat,void*),void *ctx)
2081: {
2083:   DM             dm;

2087:   KSPGetDM(ksp,&dm);
2088:   DMKSPSetComputeOperators(dm,func,ctx);
2089:   if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX;
2090:   return(0);
2091: }

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

2098:    Logically Collective

2100:    Input Arguments:
2101: +  ksp - the KSP context
2102: .  func - function to compute the right hand side
2103: -  ctx - optional context

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

2108: +  ksp - the KSP context
2109: .  b - right hand side of linear system
2110: -  ctx - optional user-provided context

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

2114:    Level: beginner

2116: .seealso: KSPSolve(), DMKSPSetComputeRHS(), KSPSetComputeOperators()
2117: @*/
2118: PetscErrorCode KSPSetComputeRHS(KSP ksp,PetscErrorCode (*func)(KSP,Vec,void*),void *ctx)
2119: {
2121:   DM             dm;

2125:   KSPGetDM(ksp,&dm);
2126:   DMKSPSetComputeRHS(dm,func,ctx);
2127:   return(0);
2128: }

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

2135:    Logically Collective

2137:    Input Arguments:
2138: +  ksp - the KSP context
2139: .  func - function to compute the initial guess
2140: -  ctx - optional context

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

2145: +  ksp - the KSP context
2146: .  x - solution vector
2147: -  ctx - optional user-provided context

2149:    Level: beginner

2151: .seealso: KSPSolve(), KSPSetComputeRHS(), KSPSetComputeOperators(), DMKSPSetComputeInitialGuess()
2152: @*/
2153: PetscErrorCode KSPSetComputeInitialGuess(KSP ksp,PetscErrorCode (*func)(KSP,Vec,void*),void *ctx)
2154: {
2156:   DM             dm;

2160:   KSPGetDM(ksp,&dm);
2161:   DMKSPSetComputeInitialGuess(dm,func,ctx);
2162:   return(0);
2163: }