Actual source code: itfunc.c

petsc-3.7.7 2017-09-25
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 excellent 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: {

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

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

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

140:    Not Collective

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

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

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

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

169:    Level: advanced

171: .keywords: KSP, compute, ritz, values

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

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

192:    Collective on KSP

194:    Input Parameter:
195: .  ksp - the KSP context

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

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

206:    Level: advanced

208: .keywords: KSP, setup, blocks

210: .seealso: PCSetUpOnBlocks(), KSPSetUp(), PCSetUp()
211: @*/
212: PetscErrorCode  KSPSetUpOnBlocks(KSP ksp)
213: {
215:   PCFailedReason pcreason;

219:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
220:   PCSetUpOnBlocks(ksp->pc);
221:   PCGetSetUpFailedReason(ksp->pc,&pcreason);
222:   if (pcreason) {
223:     ksp->reason = KSP_DIVERGED_PCSETUP_FAILED;
224:   }
225:   return(0);
226: }

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

233:    Collective on KSP

235:    Input Parameters:
236: +  ksp   - iterative context obtained from KSPCreate()
237: -  flag - PETSC_TRUE to reuse the current preconditioner

239:    Level: intermediate

241: .keywords: KSP, setup

243: .seealso: KSPCreate(), KSPSolve(), KSPDestroy(), PCSetReusePreconditioner()
244: @*/
245: PetscErrorCode  KSPSetReusePreconditioner(KSP ksp,PetscBool flag)
246: {

251:   PCSetReusePreconditioner(ksp->pc,flag);
252:   return(0);
253: }

257: /*@
258:    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

260:    Collective on KSP

262:    Input Parameters:
263: +  ksp   - iterative context obtained from KSPCreate()
264: -  flag - PETSC_TRUE to skip calling the PCSetFromOptions()

266:    Level: intermediate

268: .keywords: KSP, setup

270: .seealso: KSPCreate(), KSPSolve(), KSPDestroy(), PCSetReusePreconditioner()
271: @*/
272: PetscErrorCode  KSPSetSkipPCSetFromOptions(KSP ksp,PetscBool flag)
273: {
276:   ksp->skippcsetfromoptions = flag;
277:   return(0);
278: }

282: /*@
283:    KSPSetUp - Sets up the internal data structures for the
284:    later use of an iterative solver.

286:    Collective on KSP

288:    Input Parameter:
289: .  ksp   - iterative context obtained from KSPCreate()

291:    Level: developer

293: .keywords: KSP, setup

295: .seealso: KSPCreate(), KSPSolve(), KSPDestroy()
296: @*/
297: PetscErrorCode KSPSetUp(KSP ksp)
298: {
300:   Mat            A,B;
301:   Mat            mat,pmat;
302:   MatNullSpace   nullsp;
303:   PCFailedReason pcreason;
304: 

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

311:   if (!((PetscObject)ksp)->type_name) {
312:     KSPSetType(ksp,KSPGMRES);
313:   }
314:   KSPSetUpNorms_Private(ksp,&ksp->normtype,&ksp->pc_side);

316:   if (ksp->dmActive && !ksp->setupstage) {
317:     /* first time in so build matrix and vector data structures using DM */
318:     if (!ksp->vec_rhs) {DMCreateGlobalVector(ksp->dm,&ksp->vec_rhs);}
319:     if (!ksp->vec_sol) {DMCreateGlobalVector(ksp->dm,&ksp->vec_sol);}
320:     DMCreateMatrix(ksp->dm,&A);
321:     KSPSetOperators(ksp,A,A);
322:     PetscObjectDereference((PetscObject)A);
323:   }

325:   if (ksp->dmActive) {
326:     DMKSP kdm;
327:     DMGetDMKSP(ksp->dm,&kdm);

329:     if (kdm->ops->computeinitialguess && ksp->setupstage != KSP_SETUP_NEWRHS) {
330:       /* only computes initial guess the first time through */
331:       (*kdm->ops->computeinitialguess)(ksp,ksp->vec_sol,kdm->initialguessctx);
332:       KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
333:     }
334:     if (kdm->ops->computerhs) {
335:       (*kdm->ops->computerhs)(ksp,ksp->vec_rhs,kdm->rhsctx);
336:     }

338:     if (ksp->setupstage != KSP_SETUP_NEWRHS) {
339:       if (kdm->ops->computeoperators) {
340:         KSPGetOperators(ksp,&A,&B);
341:         (*kdm->ops->computeoperators)(ksp,A,B,kdm->operatorsctx);
342:         KSPSetOperators(ksp,A,B);
343:       } else SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONGSTATE,"You called KSPSetDM() but did not use DMKSPSetComputeOperators() or KSPSetDMActive(ksp,PETSC_FALSE);");
344:     }
345:   }

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

350:   switch (ksp->setupstage) {
351:   case KSP_SETUP_NEW:
352:     (*ksp->ops->setup)(ksp);
353:     break;
354:   case KSP_SETUP_NEWMATRIX: {   /* This should be replaced with a more general mechanism */
355:   } break;
356:   default: break;
357:   }

359:   PCGetOperators(ksp->pc,&mat,&pmat);
360:   /* scale the matrix if requested */
361:   if (ksp->dscale) {
362:     PetscScalar *xx;
363:     PetscInt    i,n;
364:     PetscBool   zeroflag = PETSC_FALSE;
365:     if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
366:     if (!ksp->diagonal) { /* allocate vector to hold diagonal */
367:       MatCreateVecs(pmat,&ksp->diagonal,0);
368:     }
369:     MatGetDiagonal(pmat,ksp->diagonal);
370:     VecGetLocalSize(ksp->diagonal,&n);
371:     VecGetArray(ksp->diagonal,&xx);
372:     for (i=0; i<n; i++) {
373:       if (xx[i] != 0.0) xx[i] = 1.0/PetscSqrtReal(PetscAbsScalar(xx[i]));
374:       else {
375:         xx[i]    = 1.0;
376:         zeroflag = PETSC_TRUE;
377:       }
378:     }
379:     VecRestoreArray(ksp->diagonal,&xx);
380:     if (zeroflag) {
381:       PetscInfo(ksp,"Zero detected in diagonal of matrix, using 1 at those locations\n");
382:     }
383:     MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
384:     if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
385:     ksp->dscalefix2 = PETSC_FALSE;
386:   }
387:   PetscLogEventEnd(KSP_SetUp,ksp,ksp->vec_rhs,ksp->vec_sol,0);
388:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
389:   PCSetErrorIfFailure(ksp->pc,ksp->errorifnotconverged);
390:   PCSetUp(ksp->pc);
391:   PCGetSetUpFailedReason(ksp->pc,&pcreason);
392:   if (pcreason) {
393:     ksp->reason = KSP_DIVERGED_PCSETUP_FAILED;
394:   }

396:   MatGetNullSpace(mat,&nullsp);
397:   if (nullsp) {
398:     PetscBool test = PETSC_FALSE;
399:     PetscOptionsGetBool(((PetscObject)ksp)->options,((PetscObject)ksp)->prefix,"-ksp_test_null_space",&test,NULL);
400:     if (test) {
401:       MatNullSpaceTest(nullsp,mat,NULL);
402:     }
403:   }
404:   ksp->setupstage = KSP_SETUP_NEWRHS;
405:   return(0);
406: }

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

413:    Collective on KSP

415:    Parameter:
416: +  ksp - iterative context obtained from KSPCreate()
417: -  viewer - the viewer to display the reason


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

423:    Level: beginner

425: .keywords: KSP, solve, linear system

427: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPConvergedDefault(),
428:           KSPSolveTranspose(), KSPGetIterationNumber()
429: @*/
430: PetscErrorCode KSPReasonView(KSP ksp,PetscViewer viewer)
431: {
433:   PetscBool      isAscii;

436:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
437:   if (isAscii) {
438:     PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
439:     if (ksp->reason > 0) {
440:       if (((PetscObject) ksp)->prefix) {
441:         PetscViewerASCIIPrintf(viewer,"Linear %s solve converged due to %s iterations %D\n",((PetscObject) ksp)->prefix,KSPConvergedReasons[ksp->reason],ksp->its);
442:       } else {
443:         PetscViewerASCIIPrintf(viewer,"Linear solve converged due to %s iterations %D\n",KSPConvergedReasons[ksp->reason],ksp->its);
444:       }
445:     } else {
446:       if (((PetscObject) ksp)->prefix) {
447:         PetscViewerASCIIPrintf(viewer,"Linear %s solve did not converge due to %s iterations %D\n",((PetscObject) ksp)->prefix,KSPConvergedReasons[ksp->reason],ksp->its);
448:       } else {
449:         PetscViewerASCIIPrintf(viewer,"Linear solve did not converge due to %s iterations %D\n",KSPConvergedReasons[ksp->reason],ksp->its);
450:       }
451:       if (ksp->reason == KSP_DIVERGED_PCSETUP_FAILED) {
452:         PCFailedReason reason;
453:         PCGetSetUpFailedReason(ksp->pc,&reason);
454:         PetscViewerASCIIPrintf(viewer,"               PCSETUP_FAILED due to %s \n",PCFailedReasons[reason]);
455:       }
456:     }
457:     PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
458:   }
459:   return(0);
460: }

463: #if defined(PETSC_HAVE_THREADSAFETY)
464: #define KSPReasonViewFromOptions KSPReasonViewFromOptionsUnsafe
466: #else
468: #endif
469: /*@C
470:   KSPReasonViewFromOptions - Processes command line options to determine if/how a KSPReason is to be viewed.

472:   Collective on KSP

474:   Input Parameters:
475: . ksp   - the KSP object

477:   Level: intermediate

479: @*/
480: PetscErrorCode KSPReasonViewFromOptions(KSP ksp)
481: {
482:   PetscErrorCode    ierr;
483:   PetscViewer       viewer;
484:   PetscBool         flg;
485:   PetscViewerFormat format;

488:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_converged_reason",&viewer,&format,&flg);
489:   if (flg) {
490:     PetscViewerPushFormat(viewer,format);
491:     KSPReasonView(ksp,viewer);
492:     PetscViewerPopFormat(viewer);
493:     PetscViewerDestroy(&viewer);
494:   }
495:   return(0);
496: }

498: #include <petscdraw.h>
501: /*@
502:    KSPSolve - Solves linear system.

504:    Collective on KSP

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

511:    Options Database Keys:
512: +  -ksp_compute_eigenvalues - compute preconditioned operators eigenvalues
513: .  -ksp_plot_eigenvalues - plot the computed eigenvalues in an X-window
514: .  -ksp_plot_eigencontours - plot the computed eigenvalues in an X-window with contours
515: .  -ksp_compute_eigenvalues_explicitly - compute the eigenvalues by forming the dense operator and using LAPACK
516: .  -ksp_plot_eigenvalues_explicitly - plot the explicitly computing eigenvalues
517: .  -ksp_view_mat binary - save matrix to the default binary viewer
518: .  -ksp_view_pmat binary - save matrix used to build preconditioner to the default binary viewer
519: .  -ksp_view_rhs binary - save right hand side vector to the default binary viewer
520: .  -ksp_view_solution binary - save computed solution vector to the default binary viewer
521:            (can be read later with src/ksp/examples/tutorials/ex10.c for testing solvers)
522: .  -ksp_view_mat_explicit - for matrix-free operators, computes the matrix entries and views them
523: .  -ksp_view_preconditioned_operator_explicit - computes the product of the preconditioner and matrix as an explicit matrix and views it
524: .  -ksp_converged_reason - print reason for converged or diverged, also prints number of iterations
525: .  -ksp_final_residual - print 2-norm of true linear system residual at the end of the solution process
526: -  -ksp_view - print the ksp data structure at the end of the system solution

528:    Notes:

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

532:    The operator is specified with KSPSetOperators().

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

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

542:    Understanding Convergence:
543:    The routines KSPMonitorSet(), KSPComputeEigenvalues(), and
544:    KSPComputeEigenvaluesExplicitly() provide information on additional
545:    options to monitor convergence and print eigenvalue information.

547:    Level: beginner

549: .keywords: KSP, solve, linear system

551: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPConvergedDefault(),
552:           KSPSolveTranspose(), KSPGetIterationNumber()
553: @*/
554: PetscErrorCode KSPSolve(KSP ksp,Vec b,Vec x)
555: {
556:   PetscErrorCode    ierr;
557:   PetscBool         flag1,flag2,flag3,flg = PETSC_FALSE,inXisinB=PETSC_FALSE,guess_zero;
558:   Mat               mat,pmat;
559:   MPI_Comm          comm;
560:   MatNullSpace      nullsp;
561:   Vec               btmp,vec_rhs=0;

567:   comm = PetscObjectComm((PetscObject)ksp);
568:   if (x && x == b) {
569:     if (!ksp->guess_zero) SETERRQ(comm,PETSC_ERR_ARG_INCOMP,"Cannot use x == b with nonzero initial guess");
570:     VecDuplicate(b,&x);
571:     inXisinB = PETSC_TRUE;
572:   }
573:   if (b) {
574:     PetscObjectReference((PetscObject)b);
575:     VecDestroy(&ksp->vec_rhs);
576:     ksp->vec_rhs = b;
577:   }
578:   if (x) {
579:     PetscObjectReference((PetscObject)x);
580:     VecDestroy(&ksp->vec_sol);
581:     ksp->vec_sol = x;
582:   }
583:   KSPViewFromOptions(ksp,NULL,"-ksp_view_pre");

585:   if (ksp->presolve) {
586:     (*ksp->presolve)(ksp,ksp->vec_rhs,ksp->vec_sol,ksp->prectx);
587:   }
588:   PetscLogEventBegin(KSP_Solve,ksp,ksp->vec_rhs,ksp->vec_sol,0);

590:   /* reset the residual history list if requested */
591:   if (ksp->res_hist_reset) ksp->res_hist_len = 0;
592:   ksp->transpose_solve = PETSC_FALSE;

594:   if (ksp->guess) {
595:     KSPFischerGuessFormGuess(ksp->guess,ksp->vec_rhs,ksp->vec_sol);
596:     ksp->guess_zero = PETSC_FALSE;
597:   }
598:   /* KSPSetUp() scales the matrix if needed */
599:   KSPSetUp(ksp);
600:   KSPSetUpOnBlocks(ksp);

602:   VecLocked(ksp->vec_sol,3);

604:   PCGetOperators(ksp->pc,&mat,&pmat);
605:   /* diagonal scale RHS if called for */
606:   if (ksp->dscale) {
607:     VecPointwiseMult(ksp->vec_rhs,ksp->vec_rhs,ksp->diagonal);
608:     /* second time in, but matrix was scaled back to original */
609:     if (ksp->dscalefix && ksp->dscalefix2) {
610:       Mat mat,pmat;

612:       PCGetOperators(ksp->pc,&mat,&pmat);
613:       MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
614:       if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
615:     }

617:     /*  scale initial guess */
618:     if (!ksp->guess_zero) {
619:       if (!ksp->truediagonal) {
620:         VecDuplicate(ksp->diagonal,&ksp->truediagonal);
621:         VecCopy(ksp->diagonal,ksp->truediagonal);
622:         VecReciprocal(ksp->truediagonal);
623:       }
624:       VecPointwiseMult(ksp->vec_sol,ksp->vec_sol,ksp->truediagonal);
625:     }
626:   }
627:   PCPreSolve(ksp->pc,ksp);

629:   if (ksp->guess_zero) { VecSet(ksp->vec_sol,0.0);}
630:   if (ksp->guess_knoll) {
631:     PCApply(ksp->pc,ksp->vec_rhs,ksp->vec_sol);
632:     KSP_RemoveNullSpace(ksp,ksp->vec_sol);
633:     ksp->guess_zero = PETSC_FALSE;
634:   }

636:   /* can we mark the initial guess as zero for this solve? */
637:   guess_zero = ksp->guess_zero;
638:   if (!ksp->guess_zero) {
639:     PetscReal norm;

641:     VecNormAvailable(ksp->vec_sol,NORM_2,&flg,&norm);
642:     if (flg && !norm) ksp->guess_zero = PETSC_TRUE;
643:   }
644:   MatGetTransposeNullSpace(pmat,&nullsp);
645:   if (nullsp) {
646:     VecDuplicate(ksp->vec_rhs,&btmp);
647:     VecCopy(ksp->vec_rhs,btmp);
648:     MatNullSpaceRemove(nullsp,btmp);
649:     vec_rhs      = ksp->vec_rhs;
650:     ksp->vec_rhs = btmp;
651:   }
652:   VecLockPush(ksp->vec_rhs);
653:   if (ksp->reason == KSP_DIVERGED_PCSETUP_FAILED) {
654:     VecSetInf(ksp->vec_sol);
655:   }
656:   (*ksp->ops->solve)(ksp);
657: 
658:   VecLockPop(ksp->vec_rhs);
659:   if (nullsp) {
660:     ksp->vec_rhs = vec_rhs;
661:     VecDestroy(&btmp);
662:   }

664:   ksp->guess_zero = guess_zero;


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

670:   KSPReasonViewFromOptions(ksp);
671:   PCPostSolve(ksp->pc,ksp);

673:   /* diagonal scale solution if called for */
674:   if (ksp->dscale) {
675:     VecPointwiseMult(ksp->vec_sol,ksp->vec_sol,ksp->diagonal);
676:     /* unscale right hand side and matrix */
677:     if (ksp->dscalefix) {
678:       Mat mat,pmat;

680:       VecReciprocal(ksp->diagonal);
681:       VecPointwiseMult(ksp->vec_rhs,ksp->vec_rhs,ksp->diagonal);
682:       PCGetOperators(ksp->pc,&mat,&pmat);
683:       MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
684:       if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
685:       VecReciprocal(ksp->diagonal);
686:       ksp->dscalefix2 = PETSC_TRUE;
687:     }
688:   }
689:   PetscLogEventEnd(KSP_Solve,ksp,ksp->vec_rhs,ksp->vec_sol,0);
690:   if (ksp->postsolve) {
691:     (*ksp->postsolve)(ksp,ksp->vec_rhs,ksp->vec_sol,ksp->postctx);
692:   }

694:   if (ksp->guess) {
695:     KSPFischerGuessUpdate(ksp->guess,ksp->vec_sol);
696:   }


699:   PCGetOperators(ksp->pc,&mat,&pmat);
700:   MatViewFromOptions(mat,(PetscObject)ksp,"-ksp_view_mat");
701:   MatViewFromOptions(pmat,(PetscObject)ksp,"-ksp_view_pmat");
702:   VecViewFromOptions(ksp->vec_rhs,(PetscObject)ksp,"-ksp_view_rhs");

704:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_compute_eigenvalues",NULL,NULL,&flag1);
705:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_plot_eigenvalues",NULL,NULL,&flag2);
706:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_plot_eigencontours",NULL,NULL,&flag3);
707:   if (flag1 || flag2 || flag3) {
708:     PetscInt    nits,n,i,neig;
709:     PetscReal   *r,*c;

711:     KSPGetIterationNumber(ksp,&nits);
712:     n    = nits+2;

714:     if (!nits) {
715:       PetscPrintf(comm,"Zero iterations in solver, cannot approximate any eigenvalues\n");
716:     } else {
717:       PetscMPIInt rank;
718:       MPI_Comm_rank(comm,&rank);
719:       PetscMalloc2(n,&r,n,&c);
720:       KSPComputeEigenvalues(ksp,n,r,c,&neig);
721:       if (flag1) {
722:         PetscPrintf(comm,"Iteratively computed eigenvalues\n");
723:         for (i=0; i<neig; i++) {
724:           if (c[i] >= 0.0) {
725:             PetscPrintf(comm,"%g + %gi\n",(double)r[i],(double)c[i]);
726:           } else {
727:             PetscPrintf(comm,"%g - %gi\n",(double)r[i],-(double)c[i]);
728:           }
729:         }
730:       }
731:       if (flag2 && !rank) {
732:         PetscDraw   draw;
733:         PetscDrawSP drawsp;

735:         if (!ksp->eigviewer) {
736:           PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Iteratively Computed Eigenvalues",PETSC_DECIDE,PETSC_DECIDE,400,400,&ksp->eigviewer);
737:         }
738:         PetscViewerDrawGetDraw(ksp->eigviewer,0,&draw);
739:         PetscDrawSPCreate(draw,1,&drawsp);
740:         for (i=0; i<neig; i++) {
741:           PetscDrawSPAddPoint(drawsp,r+i,c+i);
742:         }
743:         PetscDrawSPDraw(drawsp,PETSC_TRUE);
744:         PetscDrawSPSave(drawsp);
745:         PetscDrawSPDestroy(&drawsp);
746:       }
747:       if (flag3 && !rank) {
748:         KSPPlotEigenContours_Private(ksp,neig,r,c);
749:       }
750:       PetscFree2(r,c);
751:     }
752:   }

754:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_compute_singularvalues",NULL,NULL,&flag1);
755:   if (flag1) {
756:     PetscInt nits;

758:     KSPGetIterationNumber(ksp,&nits);
759:     if (!nits) {
760:       PetscPrintf(comm,"Zero iterations in solver, cannot approximate any singular values\n");
761:     } else {
762:       PetscReal emax,emin;

764:       KSPComputeExtremeSingularValues(ksp,&emax,&emin);
765:       PetscPrintf(comm,"Iteratively computed extreme singular values: max %g min %g max/min %g\n",(double)emax,(double)emin,(double)(emax/emin));
766:     }
767:   }

769:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_compute_eigenvalues_explicitly",NULL,NULL,&flag1);
770:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_plot_eigenvalues_explicitly",NULL,NULL,&flag2);
771:   if (flag1 || flag2) {
772:     PetscInt    n,i;
773:     PetscReal   *r,*c;
774:     PetscMPIInt rank;
775:     MPI_Comm_rank(comm,&rank);
776:     VecGetSize(ksp->vec_sol,&n);
777:     PetscMalloc2(n,&r,n,&c);
778:     KSPComputeEigenvaluesExplicitly(ksp,n,r,c);
779:     if (flag1) {
780:       PetscPrintf(comm,"Explicitly computed eigenvalues\n");
781:       for (i=0; i<n; i++) {
782:         if (c[i] >= 0.0) {
783:           PetscPrintf(comm,"%g + %gi\n",(double)r[i],(double)c[i]);
784:         } else {
785:           PetscPrintf(comm,"%g - %gi\n",(double)r[i],-(double)c[i]);
786:         }
787:       }
788:     }
789:     if (flag2 && !rank) {
790:       PetscDraw   draw;
791:       PetscDrawSP drawsp;

793:       if (!ksp->eigviewer) {
794:         PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Explicitly Computed Eigenvalues",PETSC_DECIDE,PETSC_DECIDE,400,400,&ksp->eigviewer);
795:       }
796:       PetscViewerDrawGetDraw(ksp->eigviewer,0,&draw);
797:       PetscDrawSPCreate(draw,1,&drawsp);
798:       PetscDrawSPReset(drawsp);
799:       for (i=0; i<n; i++) {
800:         PetscDrawSPAddPoint(drawsp,r+i,c+i);
801:       }
802:       PetscDrawSPDraw(drawsp,PETSC_TRUE);
803:       PetscDrawSPSave(drawsp);
804:       PetscDrawSPDestroy(&drawsp);
805:     }
806:     PetscFree2(r,c);
807:   }

809:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_view_mat_explicit",NULL,NULL,&flag2);
810:   if (flag2) {
811:     Mat A,B;
812:     PCGetOperators(ksp->pc,&A,NULL);
813:     MatComputeExplicitOperator(A,&B);
814:     MatViewFromOptions(B,(PetscObject)ksp,"-ksp_view_mat_explicit");
815:     MatDestroy(&B);
816:   }
817:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_view_preconditioned_operator_explicit",NULL,NULL,&flag2);
818:   if (flag2) {
819:     Mat B;
820:     KSPComputeExplicitOperator(ksp,&B);
821:     MatViewFromOptions(B,(PetscObject)ksp,"-ksp_view_preconditioned_operator_explicit");
822:     MatDestroy(&B);
823:   }
824:   KSPViewFromOptions(ksp,NULL,"-ksp_view");

826:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_final_residual",NULL,NULL,&flg);
827:   if (flg) {
828:     Mat       A;
829:     Vec       t;
830:     PetscReal norm;
831:     if (ksp->dscale && !ksp->dscalefix) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot compute final scale with -ksp_diagonal_scale except also with -ksp_diagonal_scale_fix");
832:     PCGetOperators(ksp->pc,&A,NULL);
833:     VecDuplicate(ksp->vec_rhs,&t);
834:     KSP_MatMult(ksp,A,ksp->vec_sol,t);
835:     VecAYPX(t, -1.0, ksp->vec_rhs);
836:     VecNorm(t,NORM_2,&norm);
837:     VecDestroy(&t);
838:     PetscPrintf(comm,"KSP final norm of residual %g\n",(double)norm);
839:   }
840:   VecViewFromOptions(ksp->vec_sol,(PetscObject)ksp,"-ksp_view_solution");

842:   if (inXisinB) {
843:     VecCopy(x,b);
844:     VecDestroy(&x);
845:   }
846:   PetscObjectSAWsBlock((PetscObject)ksp);
847:   if (ksp->errorifnotconverged && ksp->reason < 0) SETERRQ(comm,PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged");
848:   return(0);
849: }

853: /*@
854:    KSPSolveTranspose - Solves the transpose of a linear system.

856:    Collective on KSP

858:    Input Parameter:
859: +  ksp - iterative context obtained from KSPCreate()
860: .  b - right hand side vector
861: -  x - solution vector

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

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

867:    Level: developer

869: .keywords: KSP, solve, linear system

871: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPConvergedDefault(),
872:           KSPSolve()
873: @*/

875: PetscErrorCode  KSPSolveTranspose(KSP ksp,Vec b,Vec x)
876: {
878:   PetscBool      inXisinB=PETSC_FALSE;

884:   if (x == b) {
885:     VecDuplicate(b,&x);
886:     inXisinB = PETSC_TRUE;
887:   }
888:   PetscObjectReference((PetscObject)b);
889:   PetscObjectReference((PetscObject)x);
890:   VecDestroy(&ksp->vec_rhs);
891:   VecDestroy(&ksp->vec_sol);

893:   ksp->vec_rhs         = b;
894:   ksp->vec_sol         = x;
895:   ksp->transpose_solve = PETSC_TRUE;

897:   KSPSetUp(ksp);
898:   if (ksp->guess_zero) { VecSet(ksp->vec_sol,0.0);}
899:   (*ksp->ops->solve)(ksp);
900:   if (!ksp->reason) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
901:   KSPReasonViewFromOptions(ksp);
902:   if (inXisinB) {
903:     VecCopy(x,b);
904:     VecDestroy(&x);
905:   }
906:   if (ksp->errorifnotconverged && ksp->reason < 0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged");
907:   return(0);
908: }

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

915:    Collective on KSP

917:    Input Parameter:
918: .  ksp - iterative context obtained from KSPCreate()

920:    Level: beginner

922: .keywords: KSP, destroy

924: .seealso: KSPCreate(), KSPSetUp(), KSPSolve()
925: @*/
926: PetscErrorCode  KSPReset(KSP ksp)
927: {

932:   if (!ksp) return(0);
933:   if (ksp->ops->reset) {
934:     (*ksp->ops->reset)(ksp);
935:   }
936:   if (ksp->pc) {PCReset(ksp->pc);}
937:   KSPFischerGuessDestroy(&ksp->guess);
938:   VecDestroyVecs(ksp->nwork,&ksp->work);
939:   VecDestroy(&ksp->vec_rhs);
940:   VecDestroy(&ksp->vec_sol);
941:   VecDestroy(&ksp->diagonal);
942:   VecDestroy(&ksp->truediagonal);

944:   ksp->setupstage = KSP_SETUP_NEW;
945:   return(0);
946: }

950: /*@
951:    KSPDestroy - Destroys KSP context.

953:    Collective on KSP

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

958:    Level: beginner

960: .keywords: KSP, destroy

962: .seealso: KSPCreate(), KSPSetUp(), KSPSolve()
963: @*/
964: PetscErrorCode  KSPDestroy(KSP *ksp)
965: {
967:   PC             pc;

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

974:   PetscObjectSAWsViewOff((PetscObject)*ksp);
975:   /*
976:    Avoid a cascading call to PCReset(ksp->pc) from the following call:
977:    PCReset() shouldn't be called from KSPDestroy() as it is unprotected by pc's
978:    refcount (and may be shared, e.g., by other ksps).
979:    */
980:   pc         = (*ksp)->pc;
981:   (*ksp)->pc = NULL;
982:   KSPReset((*ksp));
983:   (*ksp)->pc = pc;
984:     if ((*ksp)->ops->destroy) {(*(*ksp)->ops->destroy)(*ksp);}

986:   DMDestroy(&(*ksp)->dm);
987:   PCDestroy(&(*ksp)->pc);
988:   PetscFree((*ksp)->res_hist_alloc);
989:   if ((*ksp)->convergeddestroy) {
990:     (*(*ksp)->convergeddestroy)((*ksp)->cnvP);
991:   }
992:   KSPMonitorCancel((*ksp));
993:   PetscViewerDestroy(&(*ksp)->eigviewer);
994:   PetscHeaderDestroy(ksp);
995:   return(0);
996: }

1000: /*@
1001:     KSPSetPCSide - Sets the preconditioning side.

1003:     Logically Collective on KSP

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

1008:     Output Parameter:
1009: .   side - the preconditioning side, where side is one of
1010: .vb
1011:       PC_LEFT - left preconditioning (default)
1012:       PC_RIGHT - right preconditioning
1013:       PC_SYMMETRIC - symmetric preconditioning
1014: .ve

1016:     Options Database Keys:
1017: .   -ksp_pc_side <right,left,symmetric>

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

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

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

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

1030:     Level: intermediate

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

1034: .seealso: KSPGetPCSide(), KSPSetNormType(), KSPGetNormType()
1035: @*/
1036: PetscErrorCode  KSPSetPCSide(KSP ksp,PCSide side)
1037: {
1041:   ksp->pc_side = ksp->pc_side_set = side;
1042:   return(0);
1043: }

1047: /*@
1048:     KSPGetPCSide - Gets the preconditioning side.

1050:     Not Collective

1052:     Input Parameter:
1053: .   ksp - iterative context obtained from KSPCreate()

1055:     Output Parameter:
1056: .   side - the preconditioning side, where side is one of
1057: .vb
1058:       PC_LEFT - left preconditioning (default)
1059:       PC_RIGHT - right preconditioning
1060:       PC_SYMMETRIC - symmetric preconditioning
1061: .ve

1063:     Level: intermediate

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

1067: .seealso: KSPSetPCSide()
1068: @*/
1069: PetscErrorCode  KSPGetPCSide(KSP ksp,PCSide *side)
1070: {

1076:   KSPSetUpNorms_Private(ksp,&ksp->normtype,&ksp->pc_side);
1077:   *side = ksp->pc_side;
1078:   return(0);
1079: }

1083: /*@
1084:    KSPGetTolerances - Gets the relative, absolute, divergence, and maximum
1085:    iteration tolerances used by the default KSP convergence tests.

1087:    Not Collective

1089:    Input Parameter:
1090: .  ksp - the Krylov subspace context

1092:    Output Parameters:
1093: +  rtol - the relative convergence tolerance
1094: .  abstol - the absolute convergence tolerance
1095: .  dtol - the divergence tolerance
1096: -  maxits - maximum number of iterations

1098:    Notes:
1099:    The user can specify NULL for any parameter that is not needed.

1101:    Level: intermediate

1103: .keywords: KSP, get, tolerance, absolute, relative, divergence, convergence,
1104:            maximum, iterations

1106: .seealso: KSPSetTolerances()
1107: @*/
1108: PetscErrorCode  KSPGetTolerances(KSP ksp,PetscReal *rtol,PetscReal *abstol,PetscReal *dtol,PetscInt *maxits)
1109: {
1112:   if (abstol) *abstol = ksp->abstol;
1113:   if (rtol) *rtol = ksp->rtol;
1114:   if (dtol) *dtol = ksp->divtol;
1115:   if (maxits) *maxits = ksp->max_it;
1116:   return(0);
1117: }

1121: /*@
1122:    KSPSetTolerances - Sets the relative, absolute, divergence, and maximum
1123:    iteration tolerances used by the default KSP convergence testers.

1125:    Logically Collective on KSP

1127:    Input Parameters:
1128: +  ksp - the Krylov subspace context
1129: .  rtol - the relative convergence tolerance, relative decrease in the (possibly preconditioned) residual norm
1130: .  abstol - the absolute convergence tolerance   absolute size of the (possibly preconditioned) residual norm
1131: .  dtol - the divergence tolerance,   amount (possibly preconditioned) residual norm can increase before KSPConvergedDefault() concludes that the method is diverging
1132: -  maxits - maximum number of iterations to use

1134:    Options Database Keys:
1135: +  -ksp_atol <abstol> - Sets abstol
1136: .  -ksp_rtol <rtol> - Sets rtol
1137: .  -ksp_divtol <dtol> - Sets dtol
1138: -  -ksp_max_it <maxits> - Sets maxits

1140:    Notes:
1141:    Use PETSC_DEFAULT to retain the default value of any of the tolerances.

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

1146:    Level: intermediate

1148: .keywords: KSP, set, tolerance, absolute, relative, divergence,
1149:            convergence, maximum, iterations

1151: .seealso: KSPGetTolerances(), KSPConvergedDefault(), KSPSetConvergenceTest()
1152: @*/
1153: PetscErrorCode  KSPSetTolerances(KSP ksp,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt maxits)
1154: {

1162:   if (rtol != PETSC_DEFAULT) {
1163:     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);
1164:     ksp->rtol = rtol;
1165:   }
1166:   if (abstol != PETSC_DEFAULT) {
1167:     if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
1168:     ksp->abstol = abstol;
1169:   }
1170:   if (dtol != PETSC_DEFAULT) {
1171:     if (dtol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Divergence tolerance %g must be larger than 1.0",(double)dtol);
1172:     ksp->divtol = dtol;
1173:   }
1174:   if (maxits != PETSC_DEFAULT) {
1175:     if (maxits < 0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxits);
1176:     ksp->max_it = maxits;
1177:   }
1178:   return(0);
1179: }

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

1188:    Logically Collective on KSP

1190:    Input Parameters:
1191: +  ksp - iterative context obtained from KSPCreate()
1192: -  flg - PETSC_TRUE indicates the guess is non-zero, PETSC_FALSE indicates the guess is zero

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

1197:    Level: beginner

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

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

1204: .seealso: KSPGetInitialGuessNonzero(), KSPSetInitialGuessKnoll(), KSPGetInitialGuessKnoll()
1205: @*/
1206: PetscErrorCode  KSPSetInitialGuessNonzero(KSP ksp,PetscBool flg)
1207: {
1211:   ksp->guess_zero = (PetscBool) !(int)flg;
1212:   return(0);
1213: }

1217: /*@
1218:    KSPGetInitialGuessNonzero - Determines whether the KSP solver is using
1219:    a zero initial guess.

1221:    Not Collective

1223:    Input Parameter:
1224: .  ksp - iterative context obtained from KSPCreate()

1226:    Output Parameter:
1227: .  flag - PETSC_TRUE if guess is nonzero, else PETSC_FALSE

1229:    Level: intermediate

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

1233: .seealso: KSPSetInitialGuessNonzero(), KSPSetInitialGuessKnoll(), KSPGetInitialGuessKnoll()
1234: @*/
1235: PetscErrorCode  KSPGetInitialGuessNonzero(KSP ksp,PetscBool  *flag)
1236: {
1240:   if (ksp->guess_zero) *flag = PETSC_FALSE;
1241:   else *flag = PETSC_TRUE;
1242:   return(0);
1243: }

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

1250:    Logically Collective on KSP

1252:    Input Parameters:
1253: +  ksp - iterative context obtained from KSPCreate()
1254: -  flg - PETSC_TRUE indicates you want the error generated

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

1259:    Level: intermediate

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

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

1267: .seealso: KSPGetInitialGuessNonzero(), KSPSetInitialGuessKnoll(), KSPGetInitialGuessKnoll(), KSPGetErrorIfNotConverged()
1268: @*/
1269: PetscErrorCode  KSPSetErrorIfNotConverged(KSP ksp,PetscBool flg)
1270: {
1274:   ksp->errorifnotconverged = flg;
1275:   return(0);
1276: }

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

1283:    Not Collective

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

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

1291:    Level: intermediate

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

1295: .seealso: KSPSetInitialGuessNonzero(), KSPSetInitialGuessKnoll(), KSPGetInitialGuessKnoll(), KSPSetErrorIfNotConverged()
1296: @*/
1297: PetscErrorCode  KSPGetErrorIfNotConverged(KSP ksp,PetscBool  *flag)
1298: {
1302:   *flag = ksp->errorifnotconverged;
1303:   return(0);
1304: }

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

1311:    Logically Collective on KSP

1313:    Input Parameters:
1314: +  ksp - iterative context obtained from KSPCreate()
1315: -  flg - PETSC_TRUE or PETSC_FALSE

1317:    Level: advanced


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

1322: .seealso: KSPGetInitialGuessKnoll(), KSPSetInitialGuessNonzero(), KSPGetInitialGuessNonzero()
1323: @*/
1324: PetscErrorCode  KSPSetInitialGuessKnoll(KSP ksp,PetscBool flg)
1325: {
1329:   ksp->guess_knoll = flg;
1330:   return(0);
1331: }

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

1339:    Not Collective

1341:    Input Parameter:
1342: .  ksp - iterative context obtained from KSPCreate()

1344:    Output Parameter:
1345: .  flag - PETSC_TRUE if using Knoll trick, else PETSC_FALSE

1347:    Level: advanced

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

1351: .seealso: KSPSetInitialGuessKnoll(), KSPSetInitialGuessNonzero(), KSPGetInitialGuessNonzero()
1352: @*/
1353: PetscErrorCode  KSPGetInitialGuessKnoll(KSP ksp,PetscBool  *flag)
1354: {
1358:   *flag = ksp->guess_knoll;
1359:   return(0);
1360: }

1364: /*@
1365:    KSPGetComputeSingularValues - Gets the flag indicating whether the extreme singular
1366:    values will be calculated via a Lanczos or Arnoldi process as the linear
1367:    system is solved.

1369:    Not Collective

1371:    Input Parameter:
1372: .  ksp - iterative context obtained from KSPCreate()

1374:    Output Parameter:
1375: .  flg - PETSC_TRUE or PETSC_FALSE

1377:    Options Database Key:
1378: .  -ksp_monitor_singular_value - Activates KSPSetComputeSingularValues()

1380:    Notes:
1381:    Currently this option is not valid for all iterative methods.

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

1387:    Level: advanced

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

1391: .seealso: KSPComputeExtremeSingularValues(), KSPMonitorSingularValue()
1392: @*/
1393: PetscErrorCode  KSPGetComputeSingularValues(KSP ksp,PetscBool  *flg)
1394: {
1398:   *flg = ksp->calc_sings;
1399:   return(0);
1400: }

1404: /*@
1405:    KSPSetComputeSingularValues - Sets a flag so that the extreme singular
1406:    values will be calculated via a Lanczos or Arnoldi process as the linear
1407:    system is solved.

1409:    Logically Collective on KSP

1411:    Input Parameters:
1412: +  ksp - iterative context obtained from KSPCreate()
1413: -  flg - PETSC_TRUE or PETSC_FALSE

1415:    Options Database Key:
1416: .  -ksp_monitor_singular_value - Activates KSPSetComputeSingularValues()

1418:    Notes:
1419:    Currently this option is not valid for all iterative methods.

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

1425:    Level: advanced

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

1429: .seealso: KSPComputeExtremeSingularValues(), KSPMonitorSingularValue()
1430: @*/
1431: PetscErrorCode  KSPSetComputeSingularValues(KSP ksp,PetscBool flg)
1432: {
1436:   ksp->calc_sings = flg;
1437:   return(0);
1438: }

1442: /*@
1443:    KSPGetComputeEigenvalues - Gets the flag indicating that the extreme eigenvalues
1444:    values will be calculated via a Lanczos or Arnoldi process as the linear
1445:    system is solved.

1447:    Not Collective

1449:    Input Parameter:
1450: .  ksp - iterative context obtained from KSPCreate()

1452:    Output Parameter:
1453: .  flg - PETSC_TRUE or PETSC_FALSE

1455:    Notes:
1456:    Currently this option is not valid for all iterative methods.

1458:    Level: advanced

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

1462: .seealso: KSPComputeEigenvalues(), KSPComputeEigenvaluesExplicitly()
1463: @*/
1464: PetscErrorCode  KSPGetComputeEigenvalues(KSP ksp,PetscBool  *flg)
1465: {
1469:   *flg = ksp->calc_sings;
1470:   return(0);
1471: }

1475: /*@
1476:    KSPSetComputeEigenvalues - Sets a flag so that the extreme eigenvalues
1477:    values will be calculated via a Lanczos or Arnoldi process as the linear
1478:    system is solved.

1480:    Logically Collective on KSP

1482:    Input Parameters:
1483: +  ksp - iterative context obtained from KSPCreate()
1484: -  flg - PETSC_TRUE or PETSC_FALSE

1486:    Notes:
1487:    Currently this option is not valid for all iterative methods.

1489:    Level: advanced

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

1493: .seealso: KSPComputeEigenvalues(), KSPComputeEigenvaluesExplicitly()
1494: @*/
1495: PetscErrorCode  KSPSetComputeEigenvalues(KSP ksp,PetscBool flg)
1496: {
1500:   ksp->calc_sings = flg;
1501:   return(0);
1502: }

1506: /*@
1507:    KSPSetComputeRitz - Sets a flag so that the Ritz or harmonic Ritz pairs
1508:    will be calculated via a Lanczos or Arnoldi process as the linear
1509:    system is solved.

1511:    Logically Collective on KSP

1513:    Input Parameters:
1514: +  ksp - iterative context obtained from KSPCreate()
1515: -  flg - PETSC_TRUE or PETSC_FALSE

1517:    Notes:
1518:    Currently this option is only valid for the GMRES method.

1520:    Level: advanced

1522: .keywords: KSP, set, compute, ritz

1524: .seealso: KSPComputeRitz()
1525: @*/
1526: PetscErrorCode  KSPSetComputeRitz(KSP ksp, PetscBool flg)
1527: {
1531:   ksp->calc_ritz = flg;
1532:   return(0);
1533: }

1537: /*@
1538:    KSPGetRhs - Gets the right-hand-side vector for the linear system to
1539:    be solved.

1541:    Not Collective

1543:    Input Parameter:
1544: .  ksp - iterative context obtained from KSPCreate()

1546:    Output Parameter:
1547: .  r - right-hand-side vector

1549:    Level: developer

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

1553: .seealso: KSPGetSolution(), KSPSolve()
1554: @*/
1555: PetscErrorCode  KSPGetRhs(KSP ksp,Vec *r)
1556: {
1560:   *r = ksp->vec_rhs;
1561:   return(0);
1562: }

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

1571:    Not Collective

1573:    Input Parameters:
1574: .  ksp - iterative context obtained from KSPCreate()

1576:    Output Parameters:
1577: .  v - solution vector

1579:    Level: developer

1581: .keywords: KSP, get, solution

1583: .seealso: KSPGetRhs(),  KSPBuildSolution(), KSPSolve()
1584: @*/
1585: PetscErrorCode  KSPGetSolution(KSP ksp,Vec *v)
1586: {
1590:   *v = ksp->vec_sol;
1591:   return(0);
1592: }

1596: /*@
1597:    KSPSetPC - Sets the preconditioner to be used to calculate the
1598:    application of the preconditioner on a vector.

1600:    Collective on KSP

1602:    Input Parameters:
1603: +  ksp - iterative context obtained from KSPCreate()
1604: -  pc   - the preconditioner object

1606:    Notes:
1607:    Use KSPGetPC() to retrieve the preconditioner context (for example,
1608:    to free it at the end of the computations).

1610:    Level: developer

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

1614: .seealso: KSPGetPC()
1615: @*/
1616: PetscErrorCode  KSPSetPC(KSP ksp,PC pc)
1617: {

1624:   PetscObjectReference((PetscObject)pc);
1625:   PCDestroy(&ksp->pc);
1626:   ksp->pc = pc;
1627:   PetscLogObjectParent((PetscObject)ksp,(PetscObject)ksp->pc);
1628:   return(0);
1629: }

1633: /*@
1634:    KSPGetPC - Returns a pointer to the preconditioner context
1635:    set with KSPSetPC().

1637:    Not Collective

1639:    Input Parameters:
1640: .  ksp - iterative context obtained from KSPCreate()

1642:    Output Parameter:
1643: .  pc - preconditioner context

1645:    Level: developer

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

1649: .seealso: KSPSetPC()
1650: @*/
1651: PetscErrorCode  KSPGetPC(KSP ksp,PC *pc)
1652: {

1658:   if (!ksp->pc) {
1659:     PCCreate(PetscObjectComm((PetscObject)ksp),&ksp->pc);
1660:     PetscObjectIncrementTabLevel((PetscObject)ksp->pc,(PetscObject)ksp,0);
1661:     PetscLogObjectParent((PetscObject)ksp,(PetscObject)ksp->pc);
1662:   }
1663:   *pc = ksp->pc;
1664:   return(0);
1665: }

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

1672:    Collective on KSP

1674:    Input Parameters:
1675: +  ksp - iterative context obtained from KSPCreate()
1676: .  it - iteration number
1677: -  rnorm - relative norm of the residual

1679:    Notes:
1680:    This routine is called by the KSP implementations.
1681:    It does not typically need to be called by the user.

1683:    Level: developer

1685: .seealso: KSPMonitorSet()
1686: @*/
1687: PetscErrorCode KSPMonitor(KSP ksp,PetscInt it,PetscReal rnorm)
1688: {
1689:   PetscInt       i, n = ksp->numbermonitors;

1693:   for (i=0; i<n; i++) {
1694:     (*ksp->monitor[i])(ksp,it,rnorm,ksp->monitorcontext[i]);
1695:   }
1696:   return(0);
1697: }

1701: /*

1703:     Checks if two monitors are identical; if they are then it destroys the new one
1704: */
1705: PetscErrorCode PetscMonitorCompare(PetscErrorCode (*nmon)(void),void *nmctx,PetscErrorCode (*nmdestroy)(void**),PetscErrorCode (*mon)(void),void *mctx,PetscErrorCode (*mdestroy)(void**),PetscBool *identical)
1706: {
1707:   *identical = PETSC_FALSE;
1708:   if (nmon == mon && nmdestroy == mdestroy) {
1709:     if (nmctx == mctx) *identical = PETSC_TRUE;
1710:     else if (nmdestroy == (PetscErrorCode (*)(void**)) PetscViewerAndFormatDestroy) {
1711:       PetscViewerAndFormat *old = (PetscViewerAndFormat*)mctx, *newo = (PetscViewerAndFormat*)nmctx;
1712:       if (old->viewer == newo->viewer && old->format == newo->format) *identical = PETSC_TRUE;
1713:     }
1714:     if (*identical) {
1715:       if (mdestroy) {
1717:         (*mdestroy)(&nmctx);
1718:       }
1719:     }
1720:   }
1721:   return(0);
1722: }

1726: /*@C
1727:    KSPMonitorSet - Sets an ADDITIONAL function to be called at every iteration to monitor
1728:    the residual/error etc.

1730:    Logically Collective on KSP

1732:    Input Parameters:
1733: +  ksp - iterative context obtained from KSPCreate()
1734: .  monitor - pointer to function (if this is NULL, it turns off monitoring
1735: .  mctx    - [optional] context for private data for the
1736:              monitor routine (use NULL if no context is desired)
1737: -  monitordestroy - [optional] routine that frees monitor context
1738:           (may be NULL)

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

1743: +  ksp - iterative context obtained from KSPCreate()
1744: .  it - iteration number
1745: .  rnorm - (estimated) 2-norm of (preconditioned) residual
1746: -  mctx  - optional monitoring context, as set by KSPMonitorSet()

1748:    Options Database Keys:
1749: +    -ksp_monitor        - sets KSPMonitorDefault()
1750: .    -ksp_monitor_true_residual    - sets KSPMonitorTrueResidualNorm()
1751: .    -ksp_monitor_max    - sets KSPMonitorTrueResidualMaxNorm()
1752: .    -ksp_monitor_lg_residualnorm    - sets line graph monitor,
1753:                            uses KSPMonitorLGResidualNormCreate()
1754: .    -ksp_monitor_lg_true_residualnorm   - sets line graph monitor,
1755:                            uses KSPMonitorLGResidualNormCreate()
1756: .    -ksp_monitor_singular_value    - sets KSPMonitorSingularValue()
1757: -    -ksp_monitor_cancel - cancels all monitors that have
1758:                           been hardwired into a code by
1759:                           calls to KSPMonitorSet(), but
1760:                           does not cancel those set via
1761:                           the options database.

1763:    Notes:
1764:    The default is to do nothing.  To print the residual, or preconditioned
1765:    residual if KSPSetNormType(ksp,KSP_NORM_PRECONDITIONED) was called, use
1766:    KSPMonitorDefault() as the monitoring routine, with a ASCII viewer as the
1767:    context.

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

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

1775:    Level: beginner

1777: .keywords: KSP, set, monitor

1779: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(), KSPMonitorCancel()
1780: @*/
1781: PetscErrorCode  KSPMonitorSet(KSP ksp,PetscErrorCode (*monitor)(KSP,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
1782: {
1783:   PetscInt       i;
1785:   PetscBool      identical;

1789:   for (i=0; i<ksp->numbermonitors;i++) {
1790:     PetscMonitorCompare((PetscErrorCode (*)(void))monitor,mctx,monitordestroy,(PetscErrorCode (*)(void))ksp->monitor[i],ksp->monitorcontext[i],ksp->monitordestroy[i],&identical);
1791:     if (identical) return(0);
1792:   }
1793:   if (ksp->numbermonitors >= MAXKSPMONITORS) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Too many KSP monitors set");
1794:   ksp->monitor[ksp->numbermonitors]          = monitor;
1795:   ksp->monitordestroy[ksp->numbermonitors]   = monitordestroy;
1796:   ksp->monitorcontext[ksp->numbermonitors++] = (void*)mctx;
1797:   return(0);
1798: }

1802: /*@
1803:    KSPMonitorCancel - Clears all monitors for a KSP object.

1805:    Logically Collective on KSP

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

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

1815:    Level: intermediate

1817: .keywords: KSP, set, monitor

1819: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(), KSPMonitorSet()
1820: @*/
1821: PetscErrorCode  KSPMonitorCancel(KSP ksp)
1822: {
1824:   PetscInt       i;

1828:   for (i=0; i<ksp->numbermonitors; i++) {
1829:     if (ksp->monitordestroy[i]) {
1830:       (*ksp->monitordestroy[i])(&ksp->monitorcontext[i]);
1831:     }
1832:   }
1833:   ksp->numbermonitors = 0;
1834:   return(0);
1835: }

1839: /*@C
1840:    KSPGetMonitorContext - Gets the monitoring context, as set by
1841:    KSPMonitorSet() for the FIRST monitor only.

1843:    Not Collective

1845:    Input Parameter:
1846: .  ksp - iterative context obtained from KSPCreate()

1848:    Output Parameter:
1849: .  ctx - monitoring context

1851:    Level: intermediate

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

1855: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate()
1856: @*/
1857: PetscErrorCode  KSPGetMonitorContext(KSP ksp,void **ctx)
1858: {
1861:   *ctx =      (ksp->monitorcontext[0]);
1862:   return(0);
1863: }

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

1872:    Not Collective

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

1881:    Level: advanced

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

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

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

1891: .seealso: KSPGetResidualHistory()

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


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

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

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

1923:    Not Collective

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

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

1932:    Level: advanced

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

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

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

1945: .seealso: KSPGetResidualHistory()

1947: @*/
1948: PetscErrorCode  KSPGetResidualHistory(KSP ksp,PetscReal *a[],PetscInt *na)
1949: {
1952:   if (a) *a = ksp->res_hist;
1953:   if (na) *na = ksp->res_hist_len;
1954:   return(0);
1955: }

1959: /*@C
1960:    KSPSetConvergenceTest - Sets the function to be used to determine
1961:    convergence.

1963:    Logically Collective on KSP

1965:    Input Parameters:
1966: +  ksp - iterative context obtained from KSPCreate()
1967: .  converge - pointer to int function
1968: .  cctx    - context for private data for the convergence routine (may be null)
1969: -  destroy - a routine for destroying the context (may be null)

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

1974: +  ksp - iterative context obtained from KSPCreate()
1975: .  it - iteration number
1976: .  rnorm - (estimated) 2-norm of (preconditioned) residual
1977: .  reason - the reason why it has converged or diverged
1978: -  cctx  - optional convergence context, as set by KSPSetConvergenceTest()


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

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

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

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

1995:    Level: advanced

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

1999: .seealso: KSPConvergedDefault(), KSPGetConvergenceContext(), KSPSetTolerances()
2000: @*/
2001: PetscErrorCode  KSPSetConvergenceTest(KSP ksp,PetscErrorCode (*converge)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
2002: {

2007:   if (ksp->convergeddestroy) {
2008:     (*ksp->convergeddestroy)(ksp->cnvP);
2009:   }
2010:   ksp->converged        = converge;
2011:   ksp->convergeddestroy = destroy;
2012:   ksp->cnvP             = (void*)cctx;
2013:   return(0);
2014: }

2018: /*@C
2019:    KSPGetConvergenceContext - Gets the convergence context set with
2020:    KSPSetConvergenceTest().

2022:    Not Collective

2024:    Input Parameter:
2025: .  ksp - iterative context obtained from KSPCreate()

2027:    Output Parameter:
2028: .  ctx - monitoring context

2030:    Level: advanced

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

2034: .seealso: KSPConvergedDefault(), KSPSetConvergenceTest()
2035: @*/
2036: PetscErrorCode  KSPGetConvergenceContext(KSP ksp,void **ctx)
2037: {
2040:   *ctx = ksp->cnvP;
2041:   return(0);
2042: }

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

2050:    Collective on KSP

2052:    Input Parameter:
2053: .  ctx - iterative context obtained from KSPCreate()

2055:    Output Parameter:
2056:    Provide exactly one of
2057: +  v - location to stash solution.
2058: -  V - the solution is returned in this location. This vector is created
2059:        internally. This vector should NOT be destroyed by the user with
2060:        VecDestroy().

2062:    Notes:
2063:    This routine can be used in one of two ways
2064: .vb
2065:       KSPBuildSolution(ksp,NULL,&V);
2066:    or
2067:       KSPBuildSolution(ksp,v,NULL); or KSPBuildSolution(ksp,v,&v);
2068: .ve
2069:    In the first case an internal vector is allocated to store the solution
2070:    (the user cannot destroy this vector). In the second case the solution
2071:    is generated in the vector that the user provides. Note that for certain
2072:    methods, such as KSPCG, the second case requires a copy of the solution,
2073:    while in the first case the call is essentially free since it simply
2074:    returns the vector where the solution already is stored. For some methods
2075:    like GMRES this is a reasonably expensive operation and should only be
2076:    used in truly needed.

2078:    Level: advanced

2080: .keywords: KSP, build, solution

2082: .seealso: KSPGetSolution(), KSPBuildResidual()
2083: @*/
2084: PetscErrorCode  KSPBuildSolution(KSP ksp,Vec v,Vec *V)
2085: {

2090:   if (!V && !v) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONG,"Must provide either v or V");
2091:   if (!V) V = &v;
2092:   (*ksp->ops->buildsolution)(ksp,v,V);
2093:   return(0);
2094: }

2098: /*@C
2099:    KSPBuildResidual - Builds the residual in a vector provided.

2101:    Collective on KSP

2103:    Input Parameter:
2104: .  ksp - iterative context obtained from KSPCreate()

2106:    Output Parameters:
2107: +  v - optional location to stash residual.  If v is not provided,
2108:        then a location is generated.
2109: .  t - work vector.  If not provided then one is generated.
2110: -  V - the residual

2112:    Notes:
2113:    Regardless of whether or not v is provided, the residual is
2114:    returned in V.

2116:    Level: advanced

2118: .keywords: KSP, build, residual

2120: .seealso: KSPBuildSolution()
2121: @*/
2122: PetscErrorCode  KSPBuildResidual(KSP ksp,Vec t,Vec v,Vec *V)
2123: {
2125:   PetscBool      flag = PETSC_FALSE;
2126:   Vec            w    = v,tt = t;

2130:   if (!w) {
2131:     VecDuplicate(ksp->vec_rhs,&w);
2132:     PetscLogObjectParent((PetscObject)ksp,(PetscObject)w);
2133:   }
2134:   if (!tt) {
2135:     VecDuplicate(ksp->vec_sol,&tt); flag = PETSC_TRUE;
2136:     PetscLogObjectParent((PetscObject)ksp,(PetscObject)tt);
2137:   }
2138:   (*ksp->ops->buildresidual)(ksp,tt,w,V);
2139:   if (flag) {VecDestroy(&tt);}
2140:   return(0);
2141: }

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

2149:    Logically Collective on KSP

2151:    Input Parameter:
2152: +  ksp - the KSP context
2153: -  scale - PETSC_TRUE or PETSC_FALSE

2155:    Options Database Key:
2156: +   -ksp_diagonal_scale -
2157: -   -ksp_diagonal_scale_fix - scale the matrix back AFTER the solve


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

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

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

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

2174:    Level: intermediate

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

2178: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScaleFix()
2179: @*/
2180: PetscErrorCode  KSPSetDiagonalScale(KSP ksp,PetscBool scale)
2181: {
2185:   ksp->dscale = scale;
2186:   return(0);
2187: }

2191: /*@
2192:    KSPGetDiagonalScale - Checks if KSP solver scales the matrix and
2193:                           right hand side

2195:    Not Collective

2197:    Input Parameter:
2198: .  ksp - the KSP context

2200:    Output Parameter:
2201: .  scale - PETSC_TRUE or PETSC_FALSE

2203:    Notes:
2204:     BE CAREFUL with this routine: it actually scales the matrix and right
2205:     hand side that define the system. After the system is solved the matrix
2206:     and right hand side remain scaled  unless you use KSPSetDiagonalScaleFix()

2208:    Level: intermediate

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

2212: .seealso: KSPSetDiagonalScale(), KSPSetDiagonalScaleFix()
2213: @*/
2214: PetscErrorCode  KSPGetDiagonalScale(KSP ksp,PetscBool  *scale)
2215: {
2219:   *scale = ksp->dscale;
2220:   return(0);
2221: }

2225: /*@
2226:    KSPSetDiagonalScaleFix - Tells KSP to diagonally scale the system
2227:      back after solving.

2229:    Logically Collective on KSP

2231:    Input Parameter:
2232: +  ksp - the KSP context
2233: -  fix - PETSC_TRUE to scale back after the system solve, PETSC_FALSE to not
2234:          rescale (default)

2236:    Notes:
2237:      Must be called after KSPSetDiagonalScale()

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

2244:    Level: intermediate

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

2248: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScale(), KSPGetDiagonalScaleFix()
2249: @*/
2250: PetscErrorCode  KSPSetDiagonalScaleFix(KSP ksp,PetscBool fix)
2251: {
2255:   ksp->dscalefix = fix;
2256:   return(0);
2257: }

2261: /*@
2262:    KSPGetDiagonalScaleFix - Determines if KSP diagonally scales the system
2263:      back after solving.

2265:    Not Collective

2267:    Input Parameter:
2268: .  ksp - the KSP context

2270:    Output Parameter:
2271: .  fix - PETSC_TRUE to scale back after the system solve, PETSC_FALSE to not
2272:          rescale (default)

2274:    Notes:
2275:      Must be called after KSPSetDiagonalScale()

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

2282:    Level: intermediate

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

2286: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScale(), KSPSetDiagonalScaleFix()
2287: @*/
2288: PetscErrorCode  KSPGetDiagonalScaleFix(KSP ksp,PetscBool  *fix)
2289: {
2293:   *fix = ksp->dscalefix;
2294:   return(0);
2295: }

2299: /*@C
2300:    KSPSetComputeOperators - set routine to compute the linear operators

2302:    Logically Collective

2304:    Input Arguments:
2305: +  ksp - the KSP context
2306: .  func - function to compute the operators
2307: -  ctx - optional context

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

2312: +  ksp - the KSP context
2313: .  A - the linear operator
2314: .  B - preconditioning matrix
2315: -  ctx - optional user-provided context

2317:    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
2318:           unless either KSPSetComputeOperators() or KSPSetOperators() is called before that KSPSolve() is called.

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

2322:    Level: beginner

2324: .seealso: KSPSetOperators(), KSPSetComputeRHS(), DMKSPSetComputeOperators(), KSPSetComputeInitialGuess()
2325: @*/
2326: PetscErrorCode KSPSetComputeOperators(KSP ksp,PetscErrorCode (*func)(KSP,Mat,Mat,void*),void *ctx)
2327: {
2329:   DM             dm;

2333:   KSPGetDM(ksp,&dm);
2334:   DMKSPSetComputeOperators(dm,func,ctx);
2335:   if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX;
2336:   return(0);
2337: }

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

2344:    Logically Collective

2346:    Input Arguments:
2347: +  ksp - the KSP context
2348: .  func - function to compute the right hand side
2349: -  ctx - optional context

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

2354: +  ksp - the KSP context
2355: .  b - right hand side of linear system
2356: -  ctx - optional user-provided context

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

2360:    Level: beginner

2362: .seealso: KSPSolve(), DMKSPSetComputeRHS(), KSPSetComputeOperators()
2363: @*/
2364: PetscErrorCode KSPSetComputeRHS(KSP ksp,PetscErrorCode (*func)(KSP,Vec,void*),void *ctx)
2365: {
2367:   DM             dm;

2371:   KSPGetDM(ksp,&dm);
2372:   DMKSPSetComputeRHS(dm,func,ctx);
2373:   return(0);
2374: }

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

2381:    Logically Collective

2383:    Input Arguments:
2384: +  ksp - the KSP context
2385: .  func - function to compute the initial guess
2386: -  ctx - optional context

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

2391: +  ksp - the KSP context
2392: .  x - solution vector
2393: -  ctx - optional user-provided context

2395:    Level: beginner

2397: .seealso: KSPSolve(), KSPSetComputeRHS(), KSPSetComputeOperators(), DMKSPSetComputeInitialGuess()
2398: @*/
2399: PetscErrorCode KSPSetComputeInitialGuess(KSP ksp,PetscErrorCode (*func)(KSP,Vec,void*),void *ctx)
2400: {
2402:   DM             dm;

2406:   KSPGetDM(ksp,&dm);
2407:   DMKSPSetComputeInitialGuess(dm,func,ctx);
2408:   return(0);
2409: }