Actual source code: itfunc.c
petsc-3.8.0 2017-09-26
2: /*
3: Interface KSP routines that the user calls.
4: */
6: #include <petsc/private/kspimpl.h>
7: #include <petscdm.h>
9: /*@
10: KSPComputeExtremeSingularValues - Computes the extreme singular values
11: for the preconditioned operator. Called after or during KSPSolve().
13: Not Collective
15: Input Parameter:
16: . ksp - iterative context obtained from KSPCreate()
18: Output Parameters:
19: . emin, emax - extreme singular values
21: Options Database Keys:
22: . -ksp_compute_singularvalues - compute extreme singular values and print when KSPSolve completes.
24: Notes:
25: One must call KSPSetComputeSingularValues() before calling KSPSetUp()
26: (or use the option -ksp_compute_eigenvalues) in order for this routine to work correctly.
28: Many users may just want to use the monitoring routine
29: KSPMonitorSingularValue() (which can be set with option -ksp_monitor_singular_value)
30: to print the extreme singular values at each iteration of the linear solve.
32: Estimates of the smallest singular value may be very inaccurate, especially if the Krylov method has not converged.
33: The largest singular value is usually accurate to within a few percent if the method has converged, but is still not
34: intended for eigenanalysis.
36: Disable restarts if using KSPGMRES, otherwise this estimate will only be using those iterations after the last
37: restart. See KSPGMRESSetRestart() for more details.
39: Level: advanced
41: .keywords: KSP, compute, extreme, singular, values
43: .seealso: KSPSetComputeSingularValues(), KSPMonitorSingularValue(), KSPComputeEigenvalues()
44: @*/
45: PetscErrorCode KSPComputeExtremeSingularValues(KSP ksp,PetscReal *emax,PetscReal *emin)
46: {
53: if (!ksp->calc_sings) SETERRQ(PetscObjectComm((PetscObject)ksp),4,"Singular values not requested before KSPSetUp()");
55: if (ksp->ops->computeextremesingularvalues) {
56: (*ksp->ops->computeextremesingularvalues)(ksp,emax,emin);
57: } else {
58: *emin = -1.0;
59: *emax = -1.0;
60: }
61: return(0);
62: }
64: /*@
65: KSPComputeEigenvalues - Computes the extreme eigenvalues for the
66: preconditioned operator. Called after or during KSPSolve().
68: Not Collective
70: Input Parameter:
71: + ksp - iterative context obtained from KSPCreate()
72: - n - size of arrays r and c. The number of eigenvalues computed (neig) will, in
73: general, be less than this.
75: Output Parameters:
76: + r - real part of computed eigenvalues, provided by user with a dimension of at least n
77: . c - complex part of computed eigenvalues, provided by user with a dimension of at least n
78: - neig - actual number of eigenvalues computed (will be less than or equal to n)
80: Options Database Keys:
81: + -ksp_compute_eigenvalues - Prints eigenvalues to stdout
82: - -ksp_plot_eigenvalues - Plots eigenvalues in an x-window display
84: Notes:
85: The number of eigenvalues estimated depends on the size of the Krylov space
86: generated during the KSPSolve() ; for example, with
87: CG it corresponds to the number of CG iterations, for GMRES it is the number
88: of GMRES iterations SINCE the last restart. Any extra space in r[] and c[]
89: will be ignored.
91: KSPComputeEigenvalues() does not usually provide accurate estimates; it is
92: intended only for assistance in understanding the convergence of iterative
93: methods, not for eigenanalysis. For accurate computation of eigenvalues we recommend using
94: the excellent package SLEPc.
96: One must call KSPSetComputeEigenvalues() before calling KSPSetUp()
97: in order for this routine to work correctly.
99: Many users may just want to use the monitoring routine
100: KSPMonitorSingularValue() (which can be set with option -ksp_monitor_singular_value)
101: to print the singular values at each iteration of the linear solve.
103: Level: advanced
105: .keywords: KSP, compute, extreme, singular, values
107: .seealso: KSPSetComputeSingularValues(), KSPMonitorSingularValue(), KSPComputeExtremeSingularValues()
108: @*/
109: PetscErrorCode KSPComputeEigenvalues(KSP ksp,PetscInt n,PetscReal r[],PetscReal c[],PetscInt *neig)
110: {
117: if (n<0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Requested < 0 Eigenvalues");
119: if (!ksp->calc_sings) SETERRQ(PetscObjectComm((PetscObject)ksp),4,"Eigenvalues not requested before KSPSetUp()");
121: if (n && ksp->ops->computeeigenvalues) {
122: (*ksp->ops->computeeigenvalues)(ksp,n,r,c,neig);
123: } else {
124: *neig = 0;
125: }
126: return(0);
127: }
129: /*@
130: KSPComputeRitz - Computes the Ritz or harmonic Ritz pairs associated to the
131: smallest or largest in modulus, for the preconditioned operator.
132: Called after KSPSolve().
134: Not Collective
136: Input Parameter:
137: + ksp - iterative context obtained from KSPCreate()
138: . ritz - PETSC_TRUE or PETSC_FALSE for ritz pairs or harmonic Ritz pairs, respectively
139: . small - PETSC_TRUE or PETSC_FALSE for smallest or largest (harmonic) Ritz values, respectively
140: . nrit - number of (harmonic) Ritz pairs to compute
142: Output Parameters:
143: + nrit - actual number of computed (harmonic) Ritz pairs
144: . S - multidimensional vector with Ritz vectors
145: . tetar - real part of the Ritz values
146: . tetai - imaginary part of the Ritz values
148: Notes:
149: -For GMRES, the (harmonic) Ritz pairs are computed from the Hessenberg matrix obtained during
150: the last complete cycle, or obtained at the end of the solution if the method is stopped before
151: a restart. Then, the number of actual (harmonic) Ritz pairs computed is less or equal to the restart
152: parameter for GMRES if a complete cycle has been performed or less or equal to the number of GMRES
153: iterations.
154: -Moreover, for real matrices, the (harmonic) Ritz pairs are possibly complex-valued. In such a case,
155: the routine selects the complex (harmonic) Ritz value and its conjugate, and two successive columns of S
156: are equal to the real and the imaginary parts of the associated vectors.
157: -the (harmonic) Ritz pairs are given in order of increasing (harmonic) Ritz values in modulus
158: -this is currently not implemented when PETSc is built with complex numbers
160: One must call KSPSetComputeRitz() before calling KSPSetUp()
161: in order for this routine to work correctly.
163: Level: advanced
165: .keywords: KSP, compute, ritz, values
167: .seealso: KSPSetComputeRitz()
168: @*/
169: PetscErrorCode KSPComputeRitz(KSP ksp,PetscBool ritz,PetscBool small,PetscInt *nrit,Vec S[],PetscReal tetar[],PetscReal tetai[])
170: {
175: if (!ksp->calc_ritz) SETERRQ(PetscObjectComm((PetscObject)ksp),4,"Ritz pairs not requested before KSPSetUp()");
176: if (ksp->ops->computeritz) {(*ksp->ops->computeritz)(ksp,ritz,small,nrit,S,tetar,tetai);}
177: return(0);
178: }
179: /*@
180: KSPSetUpOnBlocks - Sets up the preconditioner for each block in
181: the block Jacobi, block Gauss-Seidel, and overlapping Schwarz
182: methods.
184: Collective on KSP
186: Input Parameter:
187: . ksp - the KSP context
189: Notes:
190: KSPSetUpOnBlocks() is a routine that the user can optinally call for
191: more precise profiling (via -log_view) of the setup phase for these
192: block preconditioners. If the user does not call KSPSetUpOnBlocks(),
193: it will automatically be called from within KSPSolve().
195: Calling KSPSetUpOnBlocks() is the same as calling PCSetUpOnBlocks()
196: on the PC context within the KSP context.
198: Level: advanced
200: .keywords: KSP, setup, blocks
202: .seealso: PCSetUpOnBlocks(), KSPSetUp(), PCSetUp()
203: @*/
204: PetscErrorCode KSPSetUpOnBlocks(KSP ksp)
205: {
207: PCFailedReason pcreason;
211: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
212: PCSetUpOnBlocks(ksp->pc);
213: PCGetSetUpFailedReason(ksp->pc,&pcreason);
214: if (pcreason) {
215: ksp->reason = KSP_DIVERGED_PCSETUP_FAILED;
216: }
217: return(0);
218: }
220: /*@
221: KSPSetReusePreconditioner - reuse the current preconditioner, do not construct a new one even if the operator changes
223: Collective on KSP
225: Input Parameters:
226: + ksp - iterative context obtained from KSPCreate()
227: - flag - PETSC_TRUE to reuse the current preconditioner
229: Level: intermediate
231: .keywords: KSP, setup
233: .seealso: KSPCreate(), KSPSolve(), KSPDestroy(), PCSetReusePreconditioner()
234: @*/
235: PetscErrorCode KSPSetReusePreconditioner(KSP ksp,PetscBool flag)
236: {
241: PCSetReusePreconditioner(ksp->pc,flag);
242: return(0);
243: }
245: /*@
246: 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
248: Collective on KSP
250: Input Parameters:
251: + ksp - iterative context obtained from KSPCreate()
252: - flag - PETSC_TRUE to skip calling the PCSetFromOptions()
254: Level: intermediate
256: .keywords: KSP, setup
258: .seealso: KSPCreate(), KSPSolve(), KSPDestroy(), PCSetReusePreconditioner()
259: @*/
260: PetscErrorCode KSPSetSkipPCSetFromOptions(KSP ksp,PetscBool flag)
261: {
264: ksp->skippcsetfromoptions = flag;
265: return(0);
266: }
268: /*@
269: KSPSetUp - Sets up the internal data structures for the
270: later use of an iterative solver.
272: Collective on KSP
274: Input Parameter:
275: . ksp - iterative context obtained from KSPCreate()
277: Level: developer
279: .keywords: KSP, setup
281: .seealso: KSPCreate(), KSPSolve(), KSPDestroy()
282: @*/
283: PetscErrorCode KSPSetUp(KSP ksp)
284: {
286: Mat A,B;
287: Mat mat,pmat;
288: MatNullSpace nullsp;
289: PCFailedReason pcreason;
290:
294: /* reset the convergence flag from the previous solves */
295: ksp->reason = KSP_CONVERGED_ITERATING;
297: if (!((PetscObject)ksp)->type_name) {
298: KSPSetType(ksp,KSPGMRES);
299: }
300: KSPSetUpNorms_Private(ksp,PETSC_TRUE,&ksp->normtype,&ksp->pc_side);
302: if (ksp->dmActive && !ksp->setupstage) {
303: /* first time in so build matrix and vector data structures using DM */
304: if (!ksp->vec_rhs) {DMCreateGlobalVector(ksp->dm,&ksp->vec_rhs);}
305: if (!ksp->vec_sol) {DMCreateGlobalVector(ksp->dm,&ksp->vec_sol);}
306: DMCreateMatrix(ksp->dm,&A);
307: KSPSetOperators(ksp,A,A);
308: PetscObjectDereference((PetscObject)A);
309: }
311: if (ksp->dmActive) {
312: DMKSP kdm;
313: DMGetDMKSP(ksp->dm,&kdm);
315: if (kdm->ops->computeinitialguess && ksp->setupstage != KSP_SETUP_NEWRHS) {
316: /* only computes initial guess the first time through */
317: (*kdm->ops->computeinitialguess)(ksp,ksp->vec_sol,kdm->initialguessctx);
318: KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
319: }
320: if (kdm->ops->computerhs) {
321: (*kdm->ops->computerhs)(ksp,ksp->vec_rhs,kdm->rhsctx);
322: }
324: if (ksp->setupstage != KSP_SETUP_NEWRHS) {
325: if (kdm->ops->computeoperators) {
326: KSPGetOperators(ksp,&A,&B);
327: (*kdm->ops->computeoperators)(ksp,A,B,kdm->operatorsctx);
328: } else SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONGSTATE,"You called KSPSetDM() but did not use DMKSPSetComputeOperators() or KSPSetDMActive(ksp,PETSC_FALSE);");
329: }
330: }
332: if (ksp->setupstage == KSP_SETUP_NEWRHS) return(0);
333: PetscLogEventBegin(KSP_SetUp,ksp,ksp->vec_rhs,ksp->vec_sol,0);
335: switch (ksp->setupstage) {
336: case KSP_SETUP_NEW:
337: (*ksp->ops->setup)(ksp);
338: break;
339: case KSP_SETUP_NEWMATRIX: { /* This should be replaced with a more general mechanism */
340: if (ksp->setupnewmatrix) {
341: (*ksp->ops->setup)(ksp);
342: }
343: } break;
344: default: break;
345: }
347: PCGetOperators(ksp->pc,&mat,&pmat);
348: /* scale the matrix if requested */
349: if (ksp->dscale) {
350: PetscScalar *xx;
351: PetscInt i,n;
352: PetscBool zeroflag = PETSC_FALSE;
353: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
354: if (!ksp->diagonal) { /* allocate vector to hold diagonal */
355: MatCreateVecs(pmat,&ksp->diagonal,0);
356: }
357: MatGetDiagonal(pmat,ksp->diagonal);
358: VecGetLocalSize(ksp->diagonal,&n);
359: VecGetArray(ksp->diagonal,&xx);
360: for (i=0; i<n; i++) {
361: if (xx[i] != 0.0) xx[i] = 1.0/PetscSqrtReal(PetscAbsScalar(xx[i]));
362: else {
363: xx[i] = 1.0;
364: zeroflag = PETSC_TRUE;
365: }
366: }
367: VecRestoreArray(ksp->diagonal,&xx);
368: if (zeroflag) {
369: PetscInfo(ksp,"Zero detected in diagonal of matrix, using 1 at those locations\n");
370: }
371: MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
372: if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
373: ksp->dscalefix2 = PETSC_FALSE;
374: }
375: PetscLogEventEnd(KSP_SetUp,ksp,ksp->vec_rhs,ksp->vec_sol,0);
376: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
377: PCSetErrorIfFailure(ksp->pc,ksp->errorifnotconverged);
378: PCSetUp(ksp->pc);
379: PCGetSetUpFailedReason(ksp->pc,&pcreason);
380: if (pcreason) {
381: ksp->reason = KSP_DIVERGED_PCSETUP_FAILED;
382: }
384: MatGetNullSpace(mat,&nullsp);
385: if (nullsp) {
386: PetscBool test = PETSC_FALSE;
387: PetscOptionsGetBool(((PetscObject)ksp)->options,((PetscObject)ksp)->prefix,"-ksp_test_null_space",&test,NULL);
388: if (test) {
389: MatNullSpaceTest(nullsp,mat,NULL);
390: }
391: }
392: ksp->setupstage = KSP_SETUP_NEWRHS;
393: return(0);
394: }
396: /*@
397: KSPReasonView - Displays the reason a KSP solve converged or diverged to a viewer
399: Collective on KSP
401: Parameter:
402: + ksp - iterative context obtained from KSPCreate()
403: - viewer - the viewer to display the reason
406: Options Database Keys:
407: . -ksp_converged_reason - print reason for converged or diverged, also prints number of iterations
409: Level: beginner
411: .keywords: KSP, solve, linear system
413: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPConvergedDefault(),
414: KSPSolveTranspose(), KSPGetIterationNumber()
415: @*/
416: PetscErrorCode KSPReasonView(KSP ksp,PetscViewer viewer)
417: {
419: PetscBool isAscii;
422: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
423: if (isAscii) {
424: PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
425: if (ksp->reason > 0) {
426: if (((PetscObject) ksp)->prefix) {
427: PetscViewerASCIIPrintf(viewer,"Linear %s solve converged due to %s iterations %D\n",((PetscObject) ksp)->prefix,KSPConvergedReasons[ksp->reason],ksp->its);
428: } else {
429: PetscViewerASCIIPrintf(viewer,"Linear solve converged due to %s iterations %D\n",KSPConvergedReasons[ksp->reason],ksp->its);
430: }
431: } else {
432: if (((PetscObject) ksp)->prefix) {
433: PetscViewerASCIIPrintf(viewer,"Linear %s solve did not converge due to %s iterations %D\n",((PetscObject) ksp)->prefix,KSPConvergedReasons[ksp->reason],ksp->its);
434: } else {
435: PetscViewerASCIIPrintf(viewer,"Linear solve did not converge due to %s iterations %D\n",KSPConvergedReasons[ksp->reason],ksp->its);
436: }
437: if (ksp->reason == KSP_DIVERGED_PCSETUP_FAILED) {
438: PCFailedReason reason;
439: PCGetSetUpFailedReason(ksp->pc,&reason);
440: PetscViewerASCIIPrintf(viewer," PCSETUP_FAILED due to %s \n",PCFailedReasons[reason]);
441: }
442: }
443: PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
444: }
445: return(0);
446: }
448: #if defined(PETSC_HAVE_THREADSAFETY)
449: #define KSPReasonViewFromOptions KSPReasonViewFromOptionsUnsafe
450: #else
451: #endif
452: /*@C
453: KSPReasonViewFromOptions - Processes command line options to determine if/how a KSPReason is to be viewed.
455: Collective on KSP
457: Input Parameters:
458: . ksp - the KSP object
460: Level: intermediate
462: @*/
463: PetscErrorCode KSPReasonViewFromOptions(KSP ksp)
464: {
465: PetscErrorCode ierr;
466: PetscViewer viewer;
467: PetscBool flg;
468: PetscViewerFormat format;
471: PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_converged_reason",&viewer,&format,&flg);
472: if (flg) {
473: PetscViewerPushFormat(viewer,format);
474: KSPReasonView(ksp,viewer);
475: PetscViewerPopFormat(viewer);
476: PetscViewerDestroy(&viewer);
477: }
478: return(0);
479: }
481: #include <petscdraw.h>
482: /*@C
483: KSPSolve - Solves linear system.
485: Collective on KSP
487: Parameter:
488: + ksp - iterative context obtained from KSPCreate()
489: . b - the right hand side vector
490: - x - the solution (this may be the same vector as b, then b will be overwritten with answer)
492: Options Database Keys:
493: + -ksp_compute_eigenvalues - compute preconditioned operators eigenvalues
494: . -ksp_plot_eigenvalues - plot the computed eigenvalues in an X-window
495: . -ksp_plot_eigencontours - plot the computed eigenvalues in an X-window with contours
496: . -ksp_compute_eigenvalues_explicitly - compute the eigenvalues by forming the dense operator and using LAPACK
497: . -ksp_plot_eigenvalues_explicitly - plot the explicitly computing eigenvalues
498: . -ksp_view_mat binary - save matrix to the default binary viewer
499: . -ksp_view_pmat binary - save matrix used to build preconditioner to the default binary viewer
500: . -ksp_view_rhs binary - save right hand side vector to the default binary viewer
501: . -ksp_view_solution binary - save computed solution vector to the default binary viewer
502: (can be read later with src/ksp/examples/tutorials/ex10.c for testing solvers)
503: . -ksp_view_mat_explicit - for matrix-free operators, computes the matrix entries and views them
504: . -ksp_view_preconditioned_operator_explicit - computes the product of the preconditioner and matrix as an explicit matrix and views it
505: . -ksp_converged_reason - print reason for converged or diverged, also prints number of iterations
506: . -ksp_final_residual - print 2-norm of true linear system residual at the end of the solution process
507: - -ksp_view - print the ksp data structure at the end of the system solution
509: Notes:
511: If one uses KSPSetDM() then x or b need not be passed. Use KSPGetSolution() to access the solution in this case.
513: The operator is specified with KSPSetOperators().
515: Call KSPGetConvergedReason() to determine if the solver converged or failed and
516: why. The number of iterations can be obtained from KSPGetIterationNumber().
518: If you provide a matrix that has a MatSetNullSpace() and MatSetTransposeNullSpace() this will use that information to solve singular systems
519: in the least squares sense with a norm minimizing solution.
520: $
521: $ A x = b where b = b_p + b_t where b_t is not in the range of A (and hence by the fundamental theorem of linear algebra is in the nullspace(A') see MatSetNullSpace()
522: $
523: $ KSP first removes b_t producing the linear system A x = b_p (which has multiple solutions) and solves this to find the ||x|| minimizing solution (and hence
524: $ it finds the solution x orthogonal to the nullspace(A). The algorithm is simply in each iteration of the Krylov method we remove the nullspace(A) from the search
525: $ direction thus the solution which is a linear combination of the search directions has no component in the nullspace(A).
526: $
527: $ We recommend always using GMRES for such singular systems.
528: $ If nullspace(A) = nullspace(A') (note symmetric matrices always satisfy this property) then both left and right preconditioning will work
529: $ If nullspace(A) != nullspace(A') then left preconditioning will work but right preconditioning may not work (or it may).
531: Developer Note: The reason we cannot always solve nullspace(A) != nullspace(A') systems with right preconditioning is because we need to remove at each iteration
532: the nullspace(AB) from the search direction. While we know the nullspace(A) the nullspace(AB) equals B^-1 times the nullspace(A) but except for trivial preconditioners
533: such as diagonal scaling we cannot apply the inverse of the preconditioner to a vector and thus cannot compute the nullspace(AB).
536: If using a direct method (e.g., via the KSP solver
537: KSPPREONLY and a preconditioner such as PCLU/PCILU),
538: then its=1. See KSPSetTolerances() and KSPConvergedDefault()
539: for more details.
541: Understanding Convergence:
542: The routines KSPMonitorSet(), KSPComputeEigenvalues(), and
543: KSPComputeEigenvaluesExplicitly() provide information on additional
544: options to monitor convergence and print eigenvalue information.
546: Level: beginner
548: .keywords: KSP, solve, linear system
550: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPConvergedDefault(),
551: KSPSolveTranspose(), KSPGetIterationNumber(), MatNullSpaceCreate(), MatSetNullSpace(), MatSetTransposeNullSpace()
552: @*/
553: PetscErrorCode KSPSolve(KSP ksp,Vec b,Vec x)
554: {
555: PetscErrorCode ierr;
556: PetscBool flag1,flag2,flag3,flg = PETSC_FALSE,inXisinB=PETSC_FALSE,guess_zero;
557: Mat mat,pmat;
558: MPI_Comm comm;
559: MatNullSpace nullsp;
560: Vec btmp,vec_rhs=0;
566: comm = PetscObjectComm((PetscObject)ksp);
567: if (x && x == b) {
568: if (!ksp->guess_zero) SETERRQ(comm,PETSC_ERR_ARG_INCOMP,"Cannot use x == b with nonzero initial guess");
569: VecDuplicate(b,&x);
570: inXisinB = PETSC_TRUE;
571: }
572: if (b) {
573: PetscObjectReference((PetscObject)b);
574: VecDestroy(&ksp->vec_rhs);
575: ksp->vec_rhs = b;
576: }
577: if (x) {
578: PetscObjectReference((PetscObject)x);
579: VecDestroy(&ksp->vec_sol);
580: ksp->vec_sol = x;
581: }
582: KSPViewFromOptions(ksp,NULL,"-ksp_view_pre");
584: if (ksp->presolve) {
585: (*ksp->presolve)(ksp,ksp->vec_rhs,ksp->vec_sol,ksp->prectx);
586: }
587: PetscLogEventBegin(KSP_Solve,ksp,ksp->vec_rhs,ksp->vec_sol,0);
589: /* reset the residual history list if requested */
590: if (ksp->res_hist_reset) ksp->res_hist_len = 0;
591: ksp->transpose_solve = PETSC_FALSE;
593: if (ksp->guess) {
594: PetscObjectState ostate,state;
596: KSPGuessSetUp(ksp->guess);
597: PetscObjectStateGet((PetscObject)ksp->vec_sol,&ostate);
598: KSPGuessFormGuess(ksp->guess,ksp->vec_rhs,ksp->vec_sol);
599: PetscObjectStateGet((PetscObject)ksp->vec_sol,&state);
600: if (state != ostate) {
601: ksp->guess_zero = PETSC_FALSE;
602: } else {
603: PetscInfo(ksp,"Using zero initial guess since the KSPGuess object did not change the vector\n");
604: ksp->guess_zero = PETSC_TRUE;
605: }
606: }
608: /* KSPSetUp() scales the matrix if needed */
609: KSPSetUp(ksp);
610: KSPSetUpOnBlocks(ksp);
612: VecLocked(ksp->vec_sol,3);
614: PCGetOperators(ksp->pc,&mat,&pmat);
615: /* diagonal scale RHS if called for */
616: if (ksp->dscale) {
617: VecPointwiseMult(ksp->vec_rhs,ksp->vec_rhs,ksp->diagonal);
618: /* second time in, but matrix was scaled back to original */
619: if (ksp->dscalefix && ksp->dscalefix2) {
620: Mat mat,pmat;
622: PCGetOperators(ksp->pc,&mat,&pmat);
623: MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
624: if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
625: }
627: /* scale initial guess */
628: if (!ksp->guess_zero) {
629: if (!ksp->truediagonal) {
630: VecDuplicate(ksp->diagonal,&ksp->truediagonal);
631: VecCopy(ksp->diagonal,ksp->truediagonal);
632: VecReciprocal(ksp->truediagonal);
633: }
634: VecPointwiseMult(ksp->vec_sol,ksp->vec_sol,ksp->truediagonal);
635: }
636: }
637: PCPreSolve(ksp->pc,ksp);
639: if (ksp->guess_zero) { VecSet(ksp->vec_sol,0.0);}
640: if (ksp->guess_knoll) { /* The Knoll trick is independent on the KSPGuess specified */
641: PCApply(ksp->pc,ksp->vec_rhs,ksp->vec_sol);
642: KSP_RemoveNullSpace(ksp,ksp->vec_sol);
643: ksp->guess_zero = PETSC_FALSE;
644: }
646: /* can we mark the initial guess as zero for this solve? */
647: guess_zero = ksp->guess_zero;
648: if (!ksp->guess_zero) {
649: PetscReal norm;
651: VecNormAvailable(ksp->vec_sol,NORM_2,&flg,&norm);
652: if (flg && !norm) ksp->guess_zero = PETSC_TRUE;
653: }
654: MatGetTransposeNullSpace(pmat,&nullsp);
655: if (nullsp) {
656: VecDuplicate(ksp->vec_rhs,&btmp);
657: VecCopy(ksp->vec_rhs,btmp);
658: MatNullSpaceRemove(nullsp,btmp);
659: vec_rhs = ksp->vec_rhs;
660: ksp->vec_rhs = btmp;
661: }
662: VecLockPush(ksp->vec_rhs);
663: if (ksp->reason == KSP_DIVERGED_PCSETUP_FAILED) {
664: VecSetInf(ksp->vec_sol);
665: }
666: (*ksp->ops->solve)(ksp);
667:
668: VecLockPop(ksp->vec_rhs);
669: if (nullsp) {
670: ksp->vec_rhs = vec_rhs;
671: VecDestroy(&btmp);
672: }
674: ksp->guess_zero = guess_zero;
677: if (!ksp->reason) SETERRQ(comm,PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
678: ksp->totalits += ksp->its;
680: KSPReasonViewFromOptions(ksp);
681: PCPostSolve(ksp->pc,ksp);
683: /* diagonal scale solution if called for */
684: if (ksp->dscale) {
685: VecPointwiseMult(ksp->vec_sol,ksp->vec_sol,ksp->diagonal);
686: /* unscale right hand side and matrix */
687: if (ksp->dscalefix) {
688: Mat mat,pmat;
690: VecReciprocal(ksp->diagonal);
691: VecPointwiseMult(ksp->vec_rhs,ksp->vec_rhs,ksp->diagonal);
692: PCGetOperators(ksp->pc,&mat,&pmat);
693: MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
694: if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
695: VecReciprocal(ksp->diagonal);
696: ksp->dscalefix2 = PETSC_TRUE;
697: }
698: }
699: PetscLogEventEnd(KSP_Solve,ksp,ksp->vec_rhs,ksp->vec_sol,0);
700: if (ksp->postsolve) {
701: (*ksp->postsolve)(ksp,ksp->vec_rhs,ksp->vec_sol,ksp->postctx);
702: }
703: if (ksp->guess) {
704: KSPGuessUpdate(ksp->guess,ksp->vec_rhs,ksp->vec_sol);
705: }
707: PCGetOperators(ksp->pc,&mat,&pmat);
708: MatViewFromOptions(mat,(PetscObject)ksp,"-ksp_view_mat");
709: MatViewFromOptions(pmat,(PetscObject)ksp,"-ksp_view_pmat");
710: VecViewFromOptions(ksp->vec_rhs,(PetscObject)ksp,"-ksp_view_rhs");
712: PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_compute_eigenvalues",NULL,NULL,&flag1);
713: PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_plot_eigenvalues",NULL,NULL,&flag2);
714: PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_plot_eigencontours",NULL,NULL,&flag3);
715: if (flag1 || flag2 || flag3) {
716: PetscInt nits,n,i,neig;
717: PetscReal *r,*c;
719: KSPGetIterationNumber(ksp,&nits);
720: n = nits+2;
722: if (!nits) {
723: PetscPrintf(comm,"Zero iterations in solver, cannot approximate any eigenvalues\n");
724: } else {
725: PetscMPIInt rank;
726: MPI_Comm_rank(comm,&rank);
727: PetscMalloc2(n,&r,n,&c);
728: KSPComputeEigenvalues(ksp,n,r,c,&neig);
729: if (flag1) {
730: PetscPrintf(comm,"Iteratively computed eigenvalues\n");
731: for (i=0; i<neig; i++) {
732: if (c[i] >= 0.0) {
733: PetscPrintf(comm,"%g + %gi\n",(double)r[i],(double)c[i]);
734: } else {
735: PetscPrintf(comm,"%g - %gi\n",(double)r[i],-(double)c[i]);
736: }
737: }
738: }
739: if (flag2 && !rank) {
740: PetscDraw draw;
741: PetscDrawSP drawsp;
743: if (!ksp->eigviewer) {
744: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Iteratively Computed Eigenvalues",PETSC_DECIDE,PETSC_DECIDE,400,400,&ksp->eigviewer);
745: }
746: PetscViewerDrawGetDraw(ksp->eigviewer,0,&draw);
747: PetscDrawSPCreate(draw,1,&drawsp);
748: for (i=0; i<neig; i++) {
749: PetscDrawSPAddPoint(drawsp,r+i,c+i);
750: }
751: PetscDrawSPDraw(drawsp,PETSC_TRUE);
752: PetscDrawSPSave(drawsp);
753: PetscDrawSPDestroy(&drawsp);
754: }
755: if (flag3 && !rank) {
756: KSPPlotEigenContours_Private(ksp,neig,r,c);
757: }
758: PetscFree2(r,c);
759: }
760: }
762: PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_compute_singularvalues",NULL,NULL,&flag1);
763: if (flag1) {
764: PetscInt nits;
766: KSPGetIterationNumber(ksp,&nits);
767: if (!nits) {
768: PetscPrintf(comm,"Zero iterations in solver, cannot approximate any singular values\n");
769: } else {
770: PetscReal emax,emin;
772: KSPComputeExtremeSingularValues(ksp,&emax,&emin);
773: PetscPrintf(comm,"Iteratively computed extreme singular values: max %g min %g max/min %g\n",(double)emax,(double)emin,(double)(emax/emin));
774: }
775: }
777: PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_compute_eigenvalues_explicitly",NULL,NULL,&flag1);
778: PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_plot_eigenvalues_explicitly",NULL,NULL,&flag2);
779: if (flag1 || flag2) {
780: PetscInt n,i;
781: PetscReal *r,*c;
782: PetscMPIInt rank;
783: MPI_Comm_rank(comm,&rank);
784: VecGetSize(ksp->vec_sol,&n);
785: PetscMalloc2(n,&r,n,&c);
786: KSPComputeEigenvaluesExplicitly(ksp,n,r,c);
787: if (flag1) {
788: PetscPrintf(comm,"Explicitly computed eigenvalues\n");
789: for (i=0; i<n; i++) {
790: if (c[i] >= 0.0) {
791: PetscPrintf(comm,"%g + %gi\n",(double)r[i],(double)c[i]);
792: } else {
793: PetscPrintf(comm,"%g - %gi\n",(double)r[i],-(double)c[i]);
794: }
795: }
796: }
797: if (flag2 && !rank) {
798: PetscDraw draw;
799: PetscDrawSP drawsp;
801: if (!ksp->eigviewer) {
802: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Explicitly Computed Eigenvalues",PETSC_DECIDE,PETSC_DECIDE,400,400,&ksp->eigviewer);
803: }
804: PetscViewerDrawGetDraw(ksp->eigviewer,0,&draw);
805: PetscDrawSPCreate(draw,1,&drawsp);
806: PetscDrawSPReset(drawsp);
807: for (i=0; i<n; i++) {
808: PetscDrawSPAddPoint(drawsp,r+i,c+i);
809: }
810: PetscDrawSPDraw(drawsp,PETSC_TRUE);
811: PetscDrawSPSave(drawsp);
812: PetscDrawSPDestroy(&drawsp);
813: }
814: PetscFree2(r,c);
815: }
817: PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_view_mat_explicit",NULL,NULL,&flag2);
818: if (flag2) {
819: Mat A,B;
820: PCGetOperators(ksp->pc,&A,NULL);
821: MatComputeExplicitOperator(A,&B);
822: MatViewFromOptions(B,(PetscObject)ksp,"-ksp_view_mat_explicit");
823: MatDestroy(&B);
824: }
825: PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_view_preconditioned_operator_explicit",NULL,NULL,&flag2);
826: if (flag2) {
827: Mat B;
828: KSPComputeExplicitOperator(ksp,&B);
829: MatViewFromOptions(B,(PetscObject)ksp,"-ksp_view_preconditioned_operator_explicit");
830: MatDestroy(&B);
831: }
832: KSPViewFromOptions(ksp,NULL,"-ksp_view");
834: PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_final_residual",NULL,NULL,&flg);
835: if (flg) {
836: Mat A;
837: Vec t;
838: PetscReal norm;
839: 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");
840: PCGetOperators(ksp->pc,&A,NULL);
841: VecDuplicate(ksp->vec_rhs,&t);
842: KSP_MatMult(ksp,A,ksp->vec_sol,t);
843: VecAYPX(t, -1.0, ksp->vec_rhs);
844: VecNorm(t,NORM_2,&norm);
845: VecDestroy(&t);
846: PetscPrintf(comm,"KSP final norm of residual %g\n",(double)norm);
847: }
848: VecViewFromOptions(ksp->vec_sol,(PetscObject)ksp,"-ksp_view_solution");
850: if (inXisinB) {
851: VecCopy(x,b);
852: VecDestroy(&x);
853: }
854: PetscObjectSAWsBlock((PetscObject)ksp);
855: if (ksp->errorifnotconverged && ksp->reason < 0) SETERRQ(comm,PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged");
856: return(0);
857: }
859: /*@
860: KSPSolveTranspose - Solves the transpose of a linear system.
862: Collective on KSP
864: Input Parameter:
865: + ksp - iterative context obtained from KSPCreate()
866: . b - right hand side vector
867: - x - solution vector
869: Notes: For complex numbers this solve the non-Hermitian transpose system.
871: This currently does NOT correctly use the null space of the operator and its transpose for solving singular systems.
873: Developer Notes: We need to implement a KSPSolveHermitianTranspose()
875: Level: developer
877: .keywords: KSP, solve, linear system
879: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPConvergedDefault(),
880: KSPSolve()
881: @*/
883: PetscErrorCode KSPSolveTranspose(KSP ksp,Vec b,Vec x)
884: {
886: PetscBool inXisinB=PETSC_FALSE;
887: Vec vec_rhs = 0,btmp;
888: Mat mat,pmat;
889: MatNullSpace nullsp;
895: if (x == b) {
896: VecDuplicate(b,&x);
897: inXisinB = PETSC_TRUE;
898: }
899: PetscObjectReference((PetscObject)b);
900: PetscObjectReference((PetscObject)x);
901: VecDestroy(&ksp->vec_rhs);
902: VecDestroy(&ksp->vec_sol);
904: ksp->vec_rhs = b;
905: ksp->vec_sol = x;
906: ksp->transpose_solve = PETSC_TRUE;
908: KSPSetUp(ksp);
909: KSPSetUpOnBlocks(ksp);
910: if (ksp->guess_zero) { VecSet(ksp->vec_sol,0.0);}
912: PCGetOperators(ksp->pc,&mat,&pmat);
913: MatGetNullSpace(pmat,&nullsp);
914: if (nullsp) {
915: VecDuplicate(ksp->vec_rhs,&btmp);
916: VecCopy(ksp->vec_rhs,btmp);
917: MatNullSpaceRemove(nullsp,btmp);
918: vec_rhs = ksp->vec_rhs;
919: ksp->vec_rhs = btmp;
920: }
922: (*ksp->ops->solve)(ksp);
923: if (nullsp) {
924: ksp->vec_rhs = vec_rhs;
925: VecDestroy(&btmp);
926: }
927: if (!ksp->reason) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
928: KSPReasonViewFromOptions(ksp);
929: if (inXisinB) {
930: VecCopy(x,b);
931: VecDestroy(&x);
932: }
933: if (ksp->errorifnotconverged && ksp->reason < 0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged");
934: return(0);
935: }
937: /*@
938: KSPReset - Resets a KSP context to the kspsetupcalled = 0 state and removes any allocated Vecs and Mats
940: Collective on KSP
942: Input Parameter:
943: . ksp - iterative context obtained from KSPCreate()
945: Level: beginner
947: .keywords: KSP, destroy
949: .seealso: KSPCreate(), KSPSetUp(), KSPSolve()
950: @*/
951: PetscErrorCode KSPReset(KSP ksp)
952: {
957: if (!ksp) return(0);
958: if (ksp->ops->reset) {
959: (*ksp->ops->reset)(ksp);
960: }
961: if (ksp->pc) {PCReset(ksp->pc);}
962: if (ksp->guess) {
963: KSPGuess guess = ksp->guess;
964: if (guess->ops->reset) { (*guess->ops->reset)(guess); }
965: }
966: VecDestroyVecs(ksp->nwork,&ksp->work);
967: VecDestroy(&ksp->vec_rhs);
968: VecDestroy(&ksp->vec_sol);
969: VecDestroy(&ksp->diagonal);
970: VecDestroy(&ksp->truediagonal);
972: ksp->setupstage = KSP_SETUP_NEW;
973: return(0);
974: }
976: /*@
977: KSPDestroy - Destroys KSP context.
979: Collective on KSP
981: Input Parameter:
982: . ksp - iterative context obtained from KSPCreate()
984: Level: beginner
986: .keywords: KSP, destroy
988: .seealso: KSPCreate(), KSPSetUp(), KSPSolve()
989: @*/
990: PetscErrorCode KSPDestroy(KSP *ksp)
991: {
993: PC pc;
996: if (!*ksp) return(0);
998: if (--((PetscObject)(*ksp))->refct > 0) {*ksp = 0; return(0);}
1000: PetscObjectSAWsViewOff((PetscObject)*ksp);
1001: /*
1002: Avoid a cascading call to PCReset(ksp->pc) from the following call:
1003: PCReset() shouldn't be called from KSPDestroy() as it is unprotected by pc's
1004: refcount (and may be shared, e.g., by other ksps).
1005: */
1006: pc = (*ksp)->pc;
1007: (*ksp)->pc = NULL;
1008: KSPReset((*ksp));
1009: (*ksp)->pc = pc;
1010: if ((*ksp)->ops->destroy) {(*(*ksp)->ops->destroy)(*ksp);}
1012: DMDestroy(&(*ksp)->dm);
1013: PCDestroy(&(*ksp)->pc);
1014: PetscFree((*ksp)->res_hist_alloc);
1015: if ((*ksp)->convergeddestroy) {
1016: (*(*ksp)->convergeddestroy)((*ksp)->cnvP);
1017: }
1018: KSPMonitorCancel((*ksp));
1019: PetscViewerDestroy(&(*ksp)->eigviewer);
1020: PetscHeaderDestroy(ksp);
1021: return(0);
1022: }
1024: /*@
1025: KSPSetPCSide - Sets the preconditioning side.
1027: Logically Collective on KSP
1029: Input Parameter:
1030: . ksp - iterative context obtained from KSPCreate()
1032: Output Parameter:
1033: . side - the preconditioning side, where side is one of
1034: .vb
1035: PC_LEFT - left preconditioning (default)
1036: PC_RIGHT - right preconditioning
1037: PC_SYMMETRIC - symmetric preconditioning
1038: .ve
1040: Options Database Keys:
1041: . -ksp_pc_side <right,left,symmetric>
1043: Notes:
1044: Left preconditioning is used by default for most Krylov methods except KSPFGMRES which only supports right preconditioning.
1046: For methods changing the side of the preconditioner changes the norm type that is used, see KSPSetNormType().
1048: Symmetric preconditioning is currently available only for the KSPQCG method. Note, however, that
1049: symmetric preconditioning can be emulated by using either right or left
1050: preconditioning and a pre or post processing step.
1052: Setting the PC side often affects the default norm type. See KSPSetNormType() for details.
1054: Level: intermediate
1056: .keywords: KSP, set, right, left, symmetric, side, preconditioner, flag
1058: .seealso: KSPGetPCSide(), KSPSetNormType(), KSPGetNormType()
1059: @*/
1060: PetscErrorCode KSPSetPCSide(KSP ksp,PCSide side)
1061: {
1065: ksp->pc_side = ksp->pc_side_set = side;
1066: return(0);
1067: }
1069: /*@
1070: KSPGetPCSide - Gets the preconditioning side.
1072: Not Collective
1074: Input Parameter:
1075: . ksp - iterative context obtained from KSPCreate()
1077: Output Parameter:
1078: . side - the preconditioning side, where side is one of
1079: .vb
1080: PC_LEFT - left preconditioning (default)
1081: PC_RIGHT - right preconditioning
1082: PC_SYMMETRIC - symmetric preconditioning
1083: .ve
1085: Level: intermediate
1087: .keywords: KSP, get, right, left, symmetric, side, preconditioner, flag
1089: .seealso: KSPSetPCSide()
1090: @*/
1091: PetscErrorCode KSPGetPCSide(KSP ksp,PCSide *side)
1092: {
1098: KSPSetUpNorms_Private(ksp,PETSC_TRUE,&ksp->normtype,&ksp->pc_side);
1099: *side = ksp->pc_side;
1100: return(0);
1101: }
1103: /*@
1104: KSPGetTolerances - Gets the relative, absolute, divergence, and maximum
1105: iteration tolerances used by the default KSP convergence tests.
1107: Not Collective
1109: Input Parameter:
1110: . ksp - the Krylov subspace context
1112: Output Parameters:
1113: + rtol - the relative convergence tolerance
1114: . abstol - the absolute convergence tolerance
1115: . dtol - the divergence tolerance
1116: - maxits - maximum number of iterations
1118: Notes:
1119: The user can specify NULL for any parameter that is not needed.
1121: Level: intermediate
1123: .keywords: KSP, get, tolerance, absolute, relative, divergence, convergence,
1124: maximum, iterations
1126: .seealso: KSPSetTolerances()
1127: @*/
1128: PetscErrorCode KSPGetTolerances(KSP ksp,PetscReal *rtol,PetscReal *abstol,PetscReal *dtol,PetscInt *maxits)
1129: {
1132: if (abstol) *abstol = ksp->abstol;
1133: if (rtol) *rtol = ksp->rtol;
1134: if (dtol) *dtol = ksp->divtol;
1135: if (maxits) *maxits = ksp->max_it;
1136: return(0);
1137: }
1139: /*@
1140: KSPSetTolerances - Sets the relative, absolute, divergence, and maximum
1141: iteration tolerances used by the default KSP convergence testers.
1143: Logically Collective on KSP
1145: Input Parameters:
1146: + ksp - the Krylov subspace context
1147: . rtol - the relative convergence tolerance, relative decrease in the (possibly preconditioned) residual norm
1148: . abstol - the absolute convergence tolerance absolute size of the (possibly preconditioned) residual norm
1149: . dtol - the divergence tolerance, amount (possibly preconditioned) residual norm can increase before KSPConvergedDefault() concludes that the method is diverging
1150: - maxits - maximum number of iterations to use
1152: Options Database Keys:
1153: + -ksp_atol <abstol> - Sets abstol
1154: . -ksp_rtol <rtol> - Sets rtol
1155: . -ksp_divtol <dtol> - Sets dtol
1156: - -ksp_max_it <maxits> - Sets maxits
1158: Notes:
1159: Use PETSC_DEFAULT to retain the default value of any of the tolerances.
1161: See KSPConvergedDefault() for details how these parameters are used in the default convergence test. See also KSPSetConvergenceTest()
1162: for setting user-defined stopping criteria.
1164: Level: intermediate
1166: .keywords: KSP, set, tolerance, absolute, relative, divergence,
1167: convergence, maximum, iterations
1169: .seealso: KSPGetTolerances(), KSPConvergedDefault(), KSPSetConvergenceTest()
1170: @*/
1171: PetscErrorCode KSPSetTolerances(KSP ksp,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt maxits)
1172: {
1180: if (rtol != PETSC_DEFAULT) {
1181: 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);
1182: ksp->rtol = rtol;
1183: }
1184: if (abstol != PETSC_DEFAULT) {
1185: if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
1186: ksp->abstol = abstol;
1187: }
1188: if (dtol != PETSC_DEFAULT) {
1189: if (dtol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Divergence tolerance %g must be larger than 1.0",(double)dtol);
1190: ksp->divtol = dtol;
1191: }
1192: if (maxits != PETSC_DEFAULT) {
1193: if (maxits < 0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxits);
1194: ksp->max_it = maxits;
1195: }
1196: return(0);
1197: }
1199: /*@
1200: KSPSetInitialGuessNonzero - Tells the iterative solver that the
1201: initial guess is nonzero; otherwise KSP assumes the initial guess
1202: is to be zero (and thus zeros it out before solving).
1204: Logically Collective on KSP
1206: Input Parameters:
1207: + ksp - iterative context obtained from KSPCreate()
1208: - flg - PETSC_TRUE indicates the guess is non-zero, PETSC_FALSE indicates the guess is zero
1210: Options database keys:
1211: . -ksp_initial_guess_nonzero : use nonzero initial guess; this takes an optional truth value (0/1/no/yes/true/false)
1213: Level: beginner
1215: Notes:
1216: If this is not called the X vector is zeroed in the call to KSPSolve().
1218: .keywords: KSP, set, initial guess, nonzero
1220: .seealso: KSPGetInitialGuessNonzero(), KSPSetGuessType(), KSPGuessType
1221: @*/
1222: PetscErrorCode KSPSetInitialGuessNonzero(KSP ksp,PetscBool flg)
1223: {
1227: ksp->guess_zero = (PetscBool) !(int)flg;
1228: return(0);
1229: }
1231: /*@
1232: KSPGetInitialGuessNonzero - Determines whether the KSP solver is using
1233: a zero initial guess.
1235: Not Collective
1237: Input Parameter:
1238: . ksp - iterative context obtained from KSPCreate()
1240: Output Parameter:
1241: . flag - PETSC_TRUE if guess is nonzero, else PETSC_FALSE
1243: Level: intermediate
1245: .keywords: KSP, set, initial guess, nonzero
1247: .seealso: KSPSetInitialGuessNonzero()
1248: @*/
1249: PetscErrorCode KSPGetInitialGuessNonzero(KSP ksp,PetscBool *flag)
1250: {
1254: if (ksp->guess_zero) *flag = PETSC_FALSE;
1255: else *flag = PETSC_TRUE;
1256: return(0);
1257: }
1259: /*@
1260: KSPSetErrorIfNotConverged - Causes KSPSolve() to generate an error if the solver has not converged.
1262: Logically Collective on KSP
1264: Input Parameters:
1265: + ksp - iterative context obtained from KSPCreate()
1266: - flg - PETSC_TRUE indicates you want the error generated
1268: Options database keys:
1269: . -ksp_error_if_not_converged : this takes an optional truth value (0/1/no/yes/true/false)
1271: Level: intermediate
1273: Notes:
1274: Normally PETSc continues if a linear solver fails to converge, you can call KSPGetConvergedReason() after a KSPSolve()
1275: to determine if it has converged.
1277: .keywords: KSP
1279: .seealso: KSPGetErrorIfNotConverged()
1280: @*/
1281: PetscErrorCode KSPSetErrorIfNotConverged(KSP ksp,PetscBool flg)
1282: {
1286: ksp->errorifnotconverged = flg;
1287: return(0);
1288: }
1290: /*@
1291: KSPGetErrorIfNotConverged - Will KSPSolve() generate an error if the solver does not converge?
1293: Not Collective
1295: Input Parameter:
1296: . ksp - iterative context obtained from KSPCreate()
1298: Output Parameter:
1299: . flag - PETSC_TRUE if it will generate an error, else PETSC_FALSE
1301: Level: intermediate
1303: .keywords: KSP
1305: .seealso: KSPSetErrorIfNotConverged()
1306: @*/
1307: PetscErrorCode KSPGetErrorIfNotConverged(KSP ksp,PetscBool *flag)
1308: {
1312: *flag = ksp->errorifnotconverged;
1313: return(0);
1314: }
1316: /*@
1317: KSPSetInitialGuessKnoll - Tells the iterative solver to use PCApply(pc,b,..) to compute the initial guess (The Knoll trick)
1319: Logically Collective on KSP
1321: Input Parameters:
1322: + ksp - iterative context obtained from KSPCreate()
1323: - flg - PETSC_TRUE or PETSC_FALSE
1325: Level: advanced
1327: Developer Note: the Knoll trick is not currently implemented using the KSPGuess class
1329: .keywords: KSP, set, initial guess, nonzero
1331: .seealso: KSPGetInitialGuessKnoll(), KSPSetInitialGuessNonzero(), KSPGetInitialGuessNonzero()
1332: @*/
1333: PetscErrorCode KSPSetInitialGuessKnoll(KSP ksp,PetscBool flg)
1334: {
1338: ksp->guess_knoll = flg;
1339: return(0);
1340: }
1342: /*@
1343: KSPGetInitialGuessKnoll - Determines whether the KSP solver is using the Knoll trick (using PCApply(pc,b,...) to compute
1344: the initial guess
1346: Not Collective
1348: Input Parameter:
1349: . ksp - iterative context obtained from KSPCreate()
1351: Output Parameter:
1352: . flag - PETSC_TRUE if using Knoll trick, else PETSC_FALSE
1354: Level: advanced
1356: .keywords: KSP, set, initial guess, nonzero
1358: .seealso: KSPSetInitialGuessKnoll(), KSPSetInitialGuessNonzero(), KSPGetInitialGuessNonzero()
1359: @*/
1360: PetscErrorCode KSPGetInitialGuessKnoll(KSP ksp,PetscBool *flag)
1361: {
1365: *flag = ksp->guess_knoll;
1366: return(0);
1367: }
1369: /*@
1370: KSPGetComputeSingularValues - Gets the flag indicating whether the extreme singular
1371: values will be calculated via a Lanczos or Arnoldi process as the linear
1372: system is solved.
1374: Not Collective
1376: Input Parameter:
1377: . ksp - iterative context obtained from KSPCreate()
1379: Output Parameter:
1380: . flg - PETSC_TRUE or PETSC_FALSE
1382: Options Database Key:
1383: . -ksp_monitor_singular_value - Activates KSPSetComputeSingularValues()
1385: Notes:
1386: Currently this option is not valid for all iterative methods.
1388: Many users may just want to use the monitoring routine
1389: KSPMonitorSingularValue() (which can be set with option -ksp_monitor_singular_value)
1390: to print the singular values at each iteration of the linear solve.
1392: Level: advanced
1394: .keywords: KSP, set, compute, singular values
1396: .seealso: KSPComputeExtremeSingularValues(), KSPMonitorSingularValue()
1397: @*/
1398: PetscErrorCode KSPGetComputeSingularValues(KSP ksp,PetscBool *flg)
1399: {
1403: *flg = ksp->calc_sings;
1404: return(0);
1405: }
1407: /*@
1408: KSPSetComputeSingularValues - Sets a flag so that the extreme singular
1409: values will be calculated via a Lanczos or Arnoldi process as the linear
1410: system is solved.
1412: Logically Collective on KSP
1414: Input Parameters:
1415: + ksp - iterative context obtained from KSPCreate()
1416: - flg - PETSC_TRUE or PETSC_FALSE
1418: Options Database Key:
1419: . -ksp_monitor_singular_value - Activates KSPSetComputeSingularValues()
1421: Notes:
1422: Currently this option is not valid for all iterative methods.
1424: Many users may just want to use the monitoring routine
1425: KSPMonitorSingularValue() (which can be set with option -ksp_monitor_singular_value)
1426: to print the singular values at each iteration of the linear solve.
1428: Level: advanced
1430: .keywords: KSP, set, compute, singular values
1432: .seealso: KSPComputeExtremeSingularValues(), KSPMonitorSingularValue()
1433: @*/
1434: PetscErrorCode KSPSetComputeSingularValues(KSP ksp,PetscBool flg)
1435: {
1439: ksp->calc_sings = flg;
1440: return(0);
1441: }
1443: /*@
1444: KSPGetComputeEigenvalues - Gets the flag indicating that the extreme eigenvalues
1445: values will be calculated via a Lanczos or Arnoldi process as the linear
1446: system is solved.
1448: Not Collective
1450: Input Parameter:
1451: . ksp - iterative context obtained from KSPCreate()
1453: Output Parameter:
1454: . flg - PETSC_TRUE or PETSC_FALSE
1456: Notes:
1457: Currently this option is not valid for all iterative methods.
1459: Level: advanced
1461: .keywords: KSP, set, compute, eigenvalues
1463: .seealso: KSPComputeEigenvalues(), KSPComputeEigenvaluesExplicitly()
1464: @*/
1465: PetscErrorCode KSPGetComputeEigenvalues(KSP ksp,PetscBool *flg)
1466: {
1470: *flg = ksp->calc_sings;
1471: return(0);
1472: }
1474: /*@
1475: KSPSetComputeEigenvalues - Sets a flag so that the extreme eigenvalues
1476: values will be calculated via a Lanczos or Arnoldi process as the linear
1477: system is solved.
1479: Logically Collective on KSP
1481: Input Parameters:
1482: + ksp - iterative context obtained from KSPCreate()
1483: - flg - PETSC_TRUE or PETSC_FALSE
1485: Notes:
1486: Currently this option is not valid for all iterative methods.
1488: Level: advanced
1490: .keywords: KSP, set, compute, eigenvalues
1492: .seealso: KSPComputeEigenvalues(), KSPComputeEigenvaluesExplicitly()
1493: @*/
1494: PetscErrorCode KSPSetComputeEigenvalues(KSP ksp,PetscBool flg)
1495: {
1499: ksp->calc_sings = flg;
1500: return(0);
1501: }
1503: /*@
1504: KSPSetComputeRitz - Sets a flag so that the Ritz or harmonic Ritz pairs
1505: will be calculated via a Lanczos or Arnoldi process as the linear
1506: system is solved.
1508: Logically Collective on KSP
1510: Input Parameters:
1511: + ksp - iterative context obtained from KSPCreate()
1512: - flg - PETSC_TRUE or PETSC_FALSE
1514: Notes:
1515: Currently this option is only valid for the GMRES method.
1517: Level: advanced
1519: .keywords: KSP, set, compute, ritz
1521: .seealso: KSPComputeRitz()
1522: @*/
1523: PetscErrorCode KSPSetComputeRitz(KSP ksp, PetscBool flg)
1524: {
1528: ksp->calc_ritz = flg;
1529: return(0);
1530: }
1532: /*@
1533: KSPGetRhs - Gets the right-hand-side vector for the linear system to
1534: be solved.
1536: Not Collective
1538: Input Parameter:
1539: . ksp - iterative context obtained from KSPCreate()
1541: Output Parameter:
1542: . r - right-hand-side vector
1544: Level: developer
1546: .keywords: KSP, get, right-hand-side, rhs
1548: .seealso: KSPGetSolution(), KSPSolve()
1549: @*/
1550: PetscErrorCode KSPGetRhs(KSP ksp,Vec *r)
1551: {
1555: *r = ksp->vec_rhs;
1556: return(0);
1557: }
1559: /*@
1560: KSPGetSolution - Gets the location of the solution for the
1561: linear system to be solved. Note that this may not be where the solution
1562: is stored during the iterative process; see KSPBuildSolution().
1564: Not Collective
1566: Input Parameters:
1567: . ksp - iterative context obtained from KSPCreate()
1569: Output Parameters:
1570: . v - solution vector
1572: Level: developer
1574: .keywords: KSP, get, solution
1576: .seealso: KSPGetRhs(), KSPBuildSolution(), KSPSolve()
1577: @*/
1578: PetscErrorCode KSPGetSolution(KSP ksp,Vec *v)
1579: {
1583: *v = ksp->vec_sol;
1584: return(0);
1585: }
1587: /*@
1588: KSPSetPC - Sets the preconditioner to be used to calculate the
1589: application of the preconditioner on a vector.
1591: Collective on KSP
1593: Input Parameters:
1594: + ksp - iterative context obtained from KSPCreate()
1595: - pc - the preconditioner object
1597: Notes:
1598: Use KSPGetPC() to retrieve the preconditioner context (for example,
1599: to free it at the end of the computations).
1601: Level: developer
1603: .keywords: KSP, set, precondition, Binv
1605: .seealso: KSPGetPC()
1606: @*/
1607: PetscErrorCode KSPSetPC(KSP ksp,PC pc)
1608: {
1615: PetscObjectReference((PetscObject)pc);
1616: PCDestroy(&ksp->pc);
1617: ksp->pc = pc;
1618: PetscLogObjectParent((PetscObject)ksp,(PetscObject)ksp->pc);
1619: return(0);
1620: }
1622: /*@
1623: KSPGetPC - Returns a pointer to the preconditioner context
1624: set with KSPSetPC().
1626: Not Collective
1628: Input Parameters:
1629: . ksp - iterative context obtained from KSPCreate()
1631: Output Parameter:
1632: . pc - preconditioner context
1634: Level: developer
1636: .keywords: KSP, get, preconditioner, Binv
1638: .seealso: KSPSetPC()
1639: @*/
1640: PetscErrorCode KSPGetPC(KSP ksp,PC *pc)
1641: {
1647: if (!ksp->pc) {
1648: PCCreate(PetscObjectComm((PetscObject)ksp),&ksp->pc);
1649: PetscObjectIncrementTabLevel((PetscObject)ksp->pc,(PetscObject)ksp,0);
1650: PetscLogObjectParent((PetscObject)ksp,(PetscObject)ksp->pc);
1651: }
1652: *pc = ksp->pc;
1653: return(0);
1654: }
1656: /*@
1657: KSPMonitor - runs the user provided monitor routines, if they exist
1659: Collective on KSP
1661: Input Parameters:
1662: + ksp - iterative context obtained from KSPCreate()
1663: . it - iteration number
1664: - rnorm - relative norm of the residual
1666: Notes:
1667: This routine is called by the KSP implementations.
1668: It does not typically need to be called by the user.
1670: Level: developer
1672: .seealso: KSPMonitorSet()
1673: @*/
1674: PetscErrorCode KSPMonitor(KSP ksp,PetscInt it,PetscReal rnorm)
1675: {
1676: PetscInt i, n = ksp->numbermonitors;
1680: for (i=0; i<n; i++) {
1681: (*ksp->monitor[i])(ksp,it,rnorm,ksp->monitorcontext[i]);
1682: }
1683: return(0);
1684: }
1686: /*
1688: Checks if two monitors are identical; if they are then it destroys the new one
1689: */
1690: PetscErrorCode PetscMonitorCompare(PetscErrorCode (*nmon)(void),void *nmctx,PetscErrorCode (*nmdestroy)(void**),PetscErrorCode (*mon)(void),void *mctx,PetscErrorCode (*mdestroy)(void**),PetscBool *identical)
1691: {
1692: *identical = PETSC_FALSE;
1693: if (nmon == mon && nmdestroy == mdestroy) {
1694: if (nmctx == mctx) *identical = PETSC_TRUE;
1695: else if (nmdestroy == (PetscErrorCode (*)(void**)) PetscViewerAndFormatDestroy) {
1696: PetscViewerAndFormat *old = (PetscViewerAndFormat*)mctx, *newo = (PetscViewerAndFormat*)nmctx;
1697: if (old->viewer == newo->viewer && old->format == newo->format) *identical = PETSC_TRUE;
1698: }
1699: if (*identical) {
1700: if (mdestroy) {
1702: (*mdestroy)(&nmctx);
1703: }
1704: }
1705: }
1706: return(0);
1707: }
1709: /*@C
1710: KSPMonitorSet - Sets an ADDITIONAL function to be called at every iteration to monitor
1711: the residual/error etc.
1713: Logically Collective on KSP
1715: Input Parameters:
1716: + ksp - iterative context obtained from KSPCreate()
1717: . monitor - pointer to function (if this is NULL, it turns off monitoring
1718: . mctx - [optional] context for private data for the
1719: monitor routine (use NULL if no context is desired)
1720: - monitordestroy - [optional] routine that frees monitor context
1721: (may be NULL)
1723: Calling Sequence of monitor:
1724: $ monitor (KSP ksp, int it, PetscReal rnorm, void *mctx)
1726: + ksp - iterative context obtained from KSPCreate()
1727: . it - iteration number
1728: . rnorm - (estimated) 2-norm of (preconditioned) residual
1729: - mctx - optional monitoring context, as set by KSPMonitorSet()
1731: Options Database Keys:
1732: + -ksp_monitor - sets KSPMonitorDefault()
1733: . -ksp_monitor_true_residual - sets KSPMonitorTrueResidualNorm()
1734: . -ksp_monitor_max - sets KSPMonitorTrueResidualMaxNorm()
1735: . -ksp_monitor_lg_residualnorm - sets line graph monitor,
1736: uses KSPMonitorLGResidualNormCreate()
1737: . -ksp_monitor_lg_true_residualnorm - sets line graph monitor,
1738: uses KSPMonitorLGResidualNormCreate()
1739: . -ksp_monitor_singular_value - sets KSPMonitorSingularValue()
1740: - -ksp_monitor_cancel - cancels all monitors that have
1741: been hardwired into a code by
1742: calls to KSPMonitorSet(), but
1743: does not cancel those set via
1744: the options database.
1746: Notes:
1747: The default is to do nothing. To print the residual, or preconditioned
1748: residual if KSPSetNormType(ksp,KSP_NORM_PRECONDITIONED) was called, use
1749: KSPMonitorDefault() as the monitoring routine, with a ASCII viewer as the
1750: context.
1752: Several different monitoring routines may be set by calling
1753: KSPMonitorSet() multiple times; all will be called in the
1754: order in which they were set.
1756: Fortran notes: Only a single monitor function can be set for each KSP object
1758: Level: beginner
1760: .keywords: KSP, set, monitor
1762: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(), KSPMonitorCancel()
1763: @*/
1764: PetscErrorCode KSPMonitorSet(KSP ksp,PetscErrorCode (*monitor)(KSP,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
1765: {
1766: PetscInt i;
1768: PetscBool identical;
1772: for (i=0; i<ksp->numbermonitors;i++) {
1773: PetscMonitorCompare((PetscErrorCode (*)(void))monitor,mctx,monitordestroy,(PetscErrorCode (*)(void))ksp->monitor[i],ksp->monitorcontext[i],ksp->monitordestroy[i],&identical);
1774: if (identical) return(0);
1775: }
1776: if (ksp->numbermonitors >= MAXKSPMONITORS) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Too many KSP monitors set");
1777: ksp->monitor[ksp->numbermonitors] = monitor;
1778: ksp->monitordestroy[ksp->numbermonitors] = monitordestroy;
1779: ksp->monitorcontext[ksp->numbermonitors++] = (void*)mctx;
1780: return(0);
1781: }
1783: /*@
1784: KSPMonitorCancel - Clears all monitors for a KSP object.
1786: Logically Collective on KSP
1788: Input Parameters:
1789: . ksp - iterative context obtained from KSPCreate()
1791: Options Database Key:
1792: . -ksp_monitor_cancel - Cancels all monitors that have
1793: been hardwired into a code by calls to KSPMonitorSet(),
1794: but does not cancel those set via the options database.
1796: Level: intermediate
1798: .keywords: KSP, set, monitor
1800: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(), KSPMonitorSet()
1801: @*/
1802: PetscErrorCode KSPMonitorCancel(KSP ksp)
1803: {
1805: PetscInt i;
1809: for (i=0; i<ksp->numbermonitors; i++) {
1810: if (ksp->monitordestroy[i]) {
1811: (*ksp->monitordestroy[i])(&ksp->monitorcontext[i]);
1812: }
1813: }
1814: ksp->numbermonitors = 0;
1815: return(0);
1816: }
1818: /*@C
1819: KSPGetMonitorContext - Gets the monitoring context, as set by
1820: KSPMonitorSet() for the FIRST monitor only.
1822: Not Collective
1824: Input Parameter:
1825: . ksp - iterative context obtained from KSPCreate()
1827: Output Parameter:
1828: . ctx - monitoring context
1830: Level: intermediate
1832: .keywords: KSP, get, monitor, context
1834: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate()
1835: @*/
1836: PetscErrorCode KSPGetMonitorContext(KSP ksp,void **ctx)
1837: {
1840: *ctx = (ksp->monitorcontext[0]);
1841: return(0);
1842: }
1844: /*@
1845: KSPSetResidualHistory - Sets the array used to hold the residual history.
1846: If set, this array will contain the residual norms computed at each
1847: iteration of the solver.
1849: Not Collective
1851: Input Parameters:
1852: + ksp - iterative context obtained from KSPCreate()
1853: . a - array to hold history
1854: . na - size of a
1855: - reset - PETSC_TRUE indicates the history counter is reset to zero
1856: for each new linear solve
1858: Level: advanced
1860: Notes: The array is NOT freed by PETSc so the user needs to keep track of
1861: it and destroy once the KSP object is destroyed.
1863: If 'a' is NULL then space is allocated for the history. If 'na' PETSC_DECIDE or PETSC_DEFAULT then a
1864: default array of length 10000 is allocated.
1866: .keywords: KSP, set, residual, history, norm
1868: .seealso: KSPGetResidualHistory()
1870: @*/
1871: PetscErrorCode KSPSetResidualHistory(KSP ksp,PetscReal a[],PetscInt na,PetscBool reset)
1872: {
1878: PetscFree(ksp->res_hist_alloc);
1879: if (na != PETSC_DECIDE && na != PETSC_DEFAULT && a) {
1880: ksp->res_hist = a;
1881: ksp->res_hist_max = na;
1882: } else {
1883: if (na != PETSC_DECIDE && na != PETSC_DEFAULT) ksp->res_hist_max = na;
1884: else ksp->res_hist_max = 10000; /* like default ksp->max_it */
1885: PetscCalloc1(ksp->res_hist_max,&ksp->res_hist_alloc);
1887: ksp->res_hist = ksp->res_hist_alloc;
1888: }
1889: ksp->res_hist_len = 0;
1890: ksp->res_hist_reset = reset;
1891: return(0);
1892: }
1894: /*@C
1895: KSPGetResidualHistory - Gets the array used to hold the residual history
1896: and the number of residuals it contains.
1898: Not Collective
1900: Input Parameter:
1901: . ksp - iterative context obtained from KSPCreate()
1903: Output Parameters:
1904: + a - pointer to array to hold history (or NULL)
1905: - na - number of used entries in a (or NULL)
1907: Level: advanced
1909: Notes:
1910: Can only be called after a KSPSetResidualHistory() otherwise a and na are set to zero
1912: The Fortran version of this routine has a calling sequence
1913: $ call KSPGetResidualHistory(KSP ksp, integer na, integer ierr)
1914: note that you have passed a Fortran array into KSPSetResidualHistory() and you need
1915: to access the residual values from this Fortran array you provided. Only the na (number of
1916: residual norms currently held) is set.
1918: .keywords: KSP, get, residual, history, norm
1920: .seealso: KSPGetResidualHistory()
1922: @*/
1923: PetscErrorCode KSPGetResidualHistory(KSP ksp,PetscReal *a[],PetscInt *na)
1924: {
1927: if (a) *a = ksp->res_hist;
1928: if (na) *na = ksp->res_hist_len;
1929: return(0);
1930: }
1932: /*@C
1933: KSPSetConvergenceTest - Sets the function to be used to determine
1934: convergence.
1936: Logically Collective on KSP
1938: Input Parameters:
1939: + ksp - iterative context obtained from KSPCreate()
1940: . converge - pointer to int function
1941: . cctx - context for private data for the convergence routine (may be null)
1942: - destroy - a routine for destroying the context (may be null)
1944: Calling sequence of converge:
1945: $ converge (KSP ksp, int it, PetscReal rnorm, KSPConvergedReason *reason,void *mctx)
1947: + ksp - iterative context obtained from KSPCreate()
1948: . it - iteration number
1949: . rnorm - (estimated) 2-norm of (preconditioned) residual
1950: . reason - the reason why it has converged or diverged
1951: - cctx - optional convergence context, as set by KSPSetConvergenceTest()
1954: Notes:
1955: Must be called after the KSP type has been set so put this after
1956: a call to KSPSetType(), or KSPSetFromOptions().
1958: The default convergence test, KSPConvergedDefault(), aborts if the
1959: residual grows to more than 10000 times the initial residual.
1961: The default is a combination of relative and absolute tolerances.
1962: The residual value that is tested may be an approximation; routines
1963: that need exact values should compute them.
1965: In the default PETSc convergence test, the precise values of reason
1966: are macros such as KSP_CONVERGED_RTOL, which are defined in petscksp.h.
1968: Level: advanced
1970: .keywords: KSP, set, convergence, test, context
1972: .seealso: KSPConvergedDefault(), KSPGetConvergenceContext(), KSPSetTolerances()
1973: @*/
1974: PetscErrorCode KSPSetConvergenceTest(KSP ksp,PetscErrorCode (*converge)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
1975: {
1980: if (ksp->convergeddestroy) {
1981: (*ksp->convergeddestroy)(ksp->cnvP);
1982: }
1983: ksp->converged = converge;
1984: ksp->convergeddestroy = destroy;
1985: ksp->cnvP = (void*)cctx;
1986: return(0);
1987: }
1989: /*@C
1990: KSPGetConvergenceContext - Gets the convergence context set with
1991: KSPSetConvergenceTest().
1993: Not Collective
1995: Input Parameter:
1996: . ksp - iterative context obtained from KSPCreate()
1998: Output Parameter:
1999: . ctx - monitoring context
2001: Level: advanced
2003: .keywords: KSP, get, convergence, test, context
2005: .seealso: KSPConvergedDefault(), KSPSetConvergenceTest()
2006: @*/
2007: PetscErrorCode KSPGetConvergenceContext(KSP ksp,void **ctx)
2008: {
2011: *ctx = ksp->cnvP;
2012: return(0);
2013: }
2015: /*@C
2016: KSPBuildSolution - Builds the approximate solution in a vector provided.
2017: This routine is NOT commonly needed (see KSPSolve()).
2019: Collective on KSP
2021: Input Parameter:
2022: . ctx - iterative context obtained from KSPCreate()
2024: Output Parameter:
2025: Provide exactly one of
2026: + v - location to stash solution.
2027: - V - the solution is returned in this location. This vector is created
2028: internally. This vector should NOT be destroyed by the user with
2029: VecDestroy().
2031: Notes:
2032: This routine can be used in one of two ways
2033: .vb
2034: KSPBuildSolution(ksp,NULL,&V);
2035: or
2036: KSPBuildSolution(ksp,v,NULL); or KSPBuildSolution(ksp,v,&v);
2037: .ve
2038: In the first case an internal vector is allocated to store the solution
2039: (the user cannot destroy this vector). In the second case the solution
2040: is generated in the vector that the user provides. Note that for certain
2041: methods, such as KSPCG, the second case requires a copy of the solution,
2042: while in the first case the call is essentially free since it simply
2043: returns the vector where the solution already is stored. For some methods
2044: like GMRES this is a reasonably expensive operation and should only be
2045: used in truly needed.
2047: Level: advanced
2049: .keywords: KSP, build, solution
2051: .seealso: KSPGetSolution(), KSPBuildResidual()
2052: @*/
2053: PetscErrorCode KSPBuildSolution(KSP ksp,Vec v,Vec *V)
2054: {
2059: if (!V && !v) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONG,"Must provide either v or V");
2060: if (!V) V = &v;
2061: (*ksp->ops->buildsolution)(ksp,v,V);
2062: return(0);
2063: }
2065: /*@C
2066: KSPBuildResidual - Builds the residual in a vector provided.
2068: Collective on KSP
2070: Input Parameter:
2071: . ksp - iterative context obtained from KSPCreate()
2073: Output Parameters:
2074: + v - optional location to stash residual. If v is not provided,
2075: then a location is generated.
2076: . t - work vector. If not provided then one is generated.
2077: - V - the residual
2079: Notes:
2080: Regardless of whether or not v is provided, the residual is
2081: returned in V.
2083: Level: advanced
2085: .keywords: KSP, build, residual
2087: .seealso: KSPBuildSolution()
2088: @*/
2089: PetscErrorCode KSPBuildResidual(KSP ksp,Vec t,Vec v,Vec *V)
2090: {
2092: PetscBool flag = PETSC_FALSE;
2093: Vec w = v,tt = t;
2097: if (!w) {
2098: VecDuplicate(ksp->vec_rhs,&w);
2099: PetscLogObjectParent((PetscObject)ksp,(PetscObject)w);
2100: }
2101: if (!tt) {
2102: VecDuplicate(ksp->vec_sol,&tt); flag = PETSC_TRUE;
2103: PetscLogObjectParent((PetscObject)ksp,(PetscObject)tt);
2104: }
2105: (*ksp->ops->buildresidual)(ksp,tt,w,V);
2106: if (flag) {VecDestroy(&tt);}
2107: return(0);
2108: }
2110: /*@
2111: KSPSetDiagonalScale - Tells KSP to symmetrically diagonally scale the system
2112: before solving. This actually CHANGES the matrix (and right hand side).
2114: Logically Collective on KSP
2116: Input Parameter:
2117: + ksp - the KSP context
2118: - scale - PETSC_TRUE or PETSC_FALSE
2120: Options Database Key:
2121: + -ksp_diagonal_scale -
2122: - -ksp_diagonal_scale_fix - scale the matrix back AFTER the solve
2125: Notes: Scales the matrix by D^(-1/2) A D^(-1/2) [D^(1/2) x ] = D^(-1/2) b
2126: where D_{ii} is 1/abs(A_{ii}) unless A_{ii} is zero and then it is 1.
2128: BE CAREFUL with this routine: it actually scales the matrix and right
2129: hand side that define the system. After the system is solved the matrix
2130: and right hand side remain scaled unless you use KSPSetDiagonalScaleFix()
2132: This should NOT be used within the SNES solves if you are using a line
2133: search.
2135: If you use this with the PCType Eisenstat preconditioner than you can
2136: use the PCEisenstatSetNoDiagonalScaling() option, or -pc_eisenstat_no_diagonal_scaling
2137: to save some unneeded, redundant flops.
2139: Level: intermediate
2141: .keywords: KSP, set, options, prefix, database
2143: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScaleFix()
2144: @*/
2145: PetscErrorCode KSPSetDiagonalScale(KSP ksp,PetscBool scale)
2146: {
2150: ksp->dscale = scale;
2151: return(0);
2152: }
2154: /*@
2155: KSPGetDiagonalScale - Checks if KSP solver scales the matrix and
2156: right hand side
2158: Not Collective
2160: Input Parameter:
2161: . ksp - the KSP context
2163: Output Parameter:
2164: . scale - PETSC_TRUE or PETSC_FALSE
2166: Notes:
2167: BE CAREFUL with this routine: it actually scales the matrix and right
2168: hand side that define the system. After the system is solved the matrix
2169: and right hand side remain scaled unless you use KSPSetDiagonalScaleFix()
2171: Level: intermediate
2173: .keywords: KSP, set, options, prefix, database
2175: .seealso: KSPSetDiagonalScale(), KSPSetDiagonalScaleFix()
2176: @*/
2177: PetscErrorCode KSPGetDiagonalScale(KSP ksp,PetscBool *scale)
2178: {
2182: *scale = ksp->dscale;
2183: return(0);
2184: }
2186: /*@
2187: KSPSetDiagonalScaleFix - Tells KSP to diagonally scale the system
2188: back after solving.
2190: Logically Collective on KSP
2192: Input Parameter:
2193: + ksp - the KSP context
2194: - fix - PETSC_TRUE to scale back after the system solve, PETSC_FALSE to not
2195: rescale (default)
2197: Notes:
2198: Must be called after KSPSetDiagonalScale()
2200: Using this will slow things down, because it rescales the matrix before and
2201: after each linear solve. This is intended mainly for testing to allow one
2202: to easily get back the original system to make sure the solution computed is
2203: accurate enough.
2205: Level: intermediate
2207: .keywords: KSP, set, options, prefix, database
2209: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScale(), KSPGetDiagonalScaleFix()
2210: @*/
2211: PetscErrorCode KSPSetDiagonalScaleFix(KSP ksp,PetscBool fix)
2212: {
2216: ksp->dscalefix = fix;
2217: return(0);
2218: }
2220: /*@
2221: KSPGetDiagonalScaleFix - Determines if KSP diagonally scales the system
2222: back after solving.
2224: Not Collective
2226: Input Parameter:
2227: . ksp - the KSP context
2229: Output Parameter:
2230: . fix - PETSC_TRUE to scale back after the system solve, PETSC_FALSE to not
2231: rescale (default)
2233: Notes:
2234: Must be called after KSPSetDiagonalScale()
2236: If PETSC_TRUE will slow things down, because it rescales the matrix before and
2237: after each linear solve. This is intended mainly for testing to allow one
2238: to easily get back the original system to make sure the solution computed is
2239: accurate enough.
2241: Level: intermediate
2243: .keywords: KSP, set, options, prefix, database
2245: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScale(), KSPSetDiagonalScaleFix()
2246: @*/
2247: PetscErrorCode KSPGetDiagonalScaleFix(KSP ksp,PetscBool *fix)
2248: {
2252: *fix = ksp->dscalefix;
2253: return(0);
2254: }
2256: /*@C
2257: KSPSetComputeOperators - set routine to compute the linear operators
2259: Logically Collective
2261: Input Arguments:
2262: + ksp - the KSP context
2263: . func - function to compute the operators
2264: - ctx - optional context
2266: Calling sequence of func:
2267: $ func(KSP ksp,Mat A,Mat B,void *ctx)
2269: + ksp - the KSP context
2270: . A - the linear operator
2271: . B - preconditioning matrix
2272: - ctx - optional user-provided context
2274: 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
2275: unless either KSPSetComputeOperators() or KSPSetOperators() is called before that KSPSolve() is called.
2277: To reuse the same preconditioner for the next KSPSolve() and not compute a new one based on the most recently computed matrix call KSPSetReusePreconditioner()
2279: Level: beginner
2281: .seealso: KSPSetOperators(), KSPSetComputeRHS(), DMKSPSetComputeOperators(), KSPSetComputeInitialGuess()
2282: @*/
2283: PetscErrorCode KSPSetComputeOperators(KSP ksp,PetscErrorCode (*func)(KSP,Mat,Mat,void*),void *ctx)
2284: {
2286: DM dm;
2290: KSPGetDM(ksp,&dm);
2291: DMKSPSetComputeOperators(dm,func,ctx);
2292: if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX;
2293: return(0);
2294: }
2296: /*@C
2297: KSPSetComputeRHS - set routine to compute the right hand side of the linear system
2299: Logically Collective
2301: Input Arguments:
2302: + ksp - the KSP context
2303: . func - function to compute the right hand side
2304: - ctx - optional context
2306: Calling sequence of func:
2307: $ func(KSP ksp,Vec b,void *ctx)
2309: + ksp - the KSP context
2310: . b - right hand side of linear system
2311: - ctx - optional user-provided context
2313: Notes: The routine you provide will be called EACH you call KSPSolve() to prepare the new right hand side for that solve
2315: Level: beginner
2317: .seealso: KSPSolve(), DMKSPSetComputeRHS(), KSPSetComputeOperators()
2318: @*/
2319: PetscErrorCode KSPSetComputeRHS(KSP ksp,PetscErrorCode (*func)(KSP,Vec,void*),void *ctx)
2320: {
2322: DM dm;
2326: KSPGetDM(ksp,&dm);
2327: DMKSPSetComputeRHS(dm,func,ctx);
2328: return(0);
2329: }
2331: /*@C
2332: KSPSetComputeInitialGuess - set routine to compute the initial guess of the linear system
2334: Logically Collective
2336: Input Arguments:
2337: + ksp - the KSP context
2338: . func - function to compute the initial guess
2339: - ctx - optional context
2341: Calling sequence of func:
2342: $ func(KSP ksp,Vec x,void *ctx)
2344: + ksp - the KSP context
2345: . x - solution vector
2346: - ctx - optional user-provided context
2348: Level: beginner
2350: .seealso: KSPSolve(), KSPSetComputeRHS(), KSPSetComputeOperators(), DMKSPSetComputeInitialGuess()
2351: @*/
2352: PetscErrorCode KSPSetComputeInitialGuess(KSP ksp,PetscErrorCode (*func)(KSP,Vec,void*),void *ctx)
2353: {
2355: DM dm;
2359: KSPGetDM(ksp,&dm);
2360: DMKSPSetComputeInitialGuess(dm,func,ctx);
2361: return(0);
2362: }