Actual source code: itfunc.c
petsc-3.7.0 2016-04-25
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(dm,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_compute_eigenvalues_explicitly - compute the eigenvalues by forming the dense operator and useing LAPACK
515: . -ksp_plot_eigenvalues_explicitly - plot the explicitly computing eigenvalues
516: . -ksp_view_mat binary - save matrix to the default binary viewer
517: . -ksp_view_pmat binary - save matrix used to build preconditioner to the default binary viewer
518: . -ksp_view_rhs binary - save right hand side vector to the default binary viewer
519: . -ksp_view_solution binary - save computed solution vector to the default binary viewer
520: (can be read later with src/ksp/examples/tutorials/ex10.c for testing solvers)
521: . -ksp_view_mat_explicit - for matrix-free operators, computes the matrix entries and views them
522: . -ksp_view_preconditioned_operator_explicit - computes the product of the preconditioner and matrix as an explicit matrix and views it
523: . -ksp_converged_reason - print reason for converged or diverged, also prints number of iterations
524: . -ksp_final_residual - print 2-norm of true linear system residual at the end of the solution process
525: - -ksp_view - print the ksp data structure at the end of the system solution
527: Notes:
529: If one uses KSPSetDM() then x or b need not be passed. Use KSPGetSolution() to access the solution in this case.
531: The operator is specified with KSPSetOperators().
533: Call KSPGetConvergedReason() to determine if the solver converged or failed and
534: why. The number of iterations can be obtained from KSPGetIterationNumber().
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()
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: KSPFischerGuessFormGuess(ksp->guess,ksp->vec_rhs,ksp->vec_sol);
595: ksp->guess_zero = PETSC_FALSE;
596: }
597: /* KSPSetUp() scales the matrix if needed */
598: KSPSetUp(ksp);
599: KSPSetUpOnBlocks(ksp);
601: VecLocked(ksp->vec_sol,3);
603: PCGetOperators(ksp->pc,&mat,&pmat);
604: /* diagonal scale RHS if called for */
605: if (ksp->dscale) {
606: VecPointwiseMult(ksp->vec_rhs,ksp->vec_rhs,ksp->diagonal);
607: /* second time in, but matrix was scaled back to original */
608: if (ksp->dscalefix && ksp->dscalefix2) {
609: Mat mat,pmat;
611: PCGetOperators(ksp->pc,&mat,&pmat);
612: MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
613: if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
614: }
616: /* scale initial guess */
617: if (!ksp->guess_zero) {
618: if (!ksp->truediagonal) {
619: VecDuplicate(ksp->diagonal,&ksp->truediagonal);
620: VecCopy(ksp->diagonal,ksp->truediagonal);
621: VecReciprocal(ksp->truediagonal);
622: }
623: VecPointwiseMult(ksp->vec_sol,ksp->vec_sol,ksp->truediagonal);
624: }
625: }
626: PCPreSolve(ksp->pc,ksp);
628: if (ksp->guess_zero) { VecSet(ksp->vec_sol,0.0);}
629: if (ksp->guess_knoll) {
630: PCApply(ksp->pc,ksp->vec_rhs,ksp->vec_sol);
631: KSP_RemoveNullSpace(ksp,ksp->vec_sol);
632: ksp->guess_zero = PETSC_FALSE;
633: }
635: /* can we mark the initial guess as zero for this solve? */
636: guess_zero = ksp->guess_zero;
637: if (!ksp->guess_zero) {
638: PetscReal norm;
640: VecNormAvailable(ksp->vec_sol,NORM_2,&flg,&norm);
641: if (flg && !norm) ksp->guess_zero = PETSC_TRUE;
642: }
643: MatGetTransposeNullSpace(pmat,&nullsp);
644: if (nullsp) {
645: VecDuplicate(ksp->vec_rhs,&btmp);
646: VecCopy(ksp->vec_rhs,btmp);
647: MatNullSpaceRemove(nullsp,btmp);
648: vec_rhs = ksp->vec_rhs;
649: ksp->vec_rhs = btmp;
650: }
651: VecLockPush(ksp->vec_rhs);
652: if (ksp->reason == KSP_DIVERGED_PCSETUP_FAILED) {
653: VecSetInf(ksp->vec_sol);
654: }
655: (*ksp->ops->solve)(ksp);
656:
657: VecLockPop(ksp->vec_rhs);
658: if (nullsp) {
659: ksp->vec_rhs = vec_rhs;
660: VecDestroy(&btmp);
661: }
663: ksp->guess_zero = guess_zero;
666: if (!ksp->reason) SETERRQ(comm,PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
667: ksp->totalits += ksp->its;
669: KSPReasonViewFromOptions(ksp);
670: PCPostSolve(ksp->pc,ksp);
672: /* diagonal scale solution if called for */
673: if (ksp->dscale) {
674: VecPointwiseMult(ksp->vec_sol,ksp->vec_sol,ksp->diagonal);
675: /* unscale right hand side and matrix */
676: if (ksp->dscalefix) {
677: Mat mat,pmat;
679: VecReciprocal(ksp->diagonal);
680: VecPointwiseMult(ksp->vec_rhs,ksp->vec_rhs,ksp->diagonal);
681: PCGetOperators(ksp->pc,&mat,&pmat);
682: MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
683: if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
684: VecReciprocal(ksp->diagonal);
685: ksp->dscalefix2 = PETSC_TRUE;
686: }
687: }
688: PetscLogEventEnd(KSP_Solve,ksp,ksp->vec_rhs,ksp->vec_sol,0);
689: if (ksp->postsolve) {
690: (*ksp->postsolve)(ksp,ksp->vec_rhs,ksp->vec_sol,ksp->postctx);
691: }
693: if (ksp->guess) {
694: KSPFischerGuessUpdate(ksp->guess,ksp->vec_sol);
695: }
698: PCGetOperators(ksp->pc,&mat,&pmat);
699: MatViewFromOptions(mat,(PetscObject)ksp,"-ksp_view_mat");
700: MatViewFromOptions(pmat,(PetscObject)ksp,"-ksp_view_pmat");
701: VecViewFromOptions(ksp->vec_rhs,(PetscObject)ksp,"-ksp_view_rhs");
703: PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_compute_eigenvalues",NULL,NULL,&flag1);
704: PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_plot_eigenvalues",NULL,NULL,&flag2);
705: PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_plot_eigenontours",NULL,NULL,&flag3);
706: if (flag1 || flag2 || flag3) {
707: PetscInt nits,n,i,neig;
708: PetscReal *r,*c;
710: KSPGetIterationNumber(ksp,&nits);
711: n = nits+2;
713: if (!nits) {
714: PetscPrintf(comm,"Zero iterations in solver, cannot approximate any eigenvalues\n");
715: } else {
716: PetscMPIInt rank;
717: MPI_Comm_rank(comm,&rank);
718: PetscMalloc2(n,&r,n,&c);
719: KSPComputeEigenvalues(ksp,n,r,c,&neig);
720: if (flag1) {
721: PetscPrintf(comm,"Iteratively computed eigenvalues\n");
722: for (i=0; i<neig; i++) {
723: if (c[i] >= 0.0) {
724: PetscPrintf(comm,"%g + %gi\n",(double)r[i],(double)c[i]);
725: } else {
726: PetscPrintf(comm,"%g - %gi\n",(double)r[i],-(double)c[i]);
727: }
728: }
729: }
730: if (flag2 && !rank) {
731: PetscDraw draw;
732: PetscDrawSP drawsp;
734: if (!ksp->eigviewer) {
735: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Iteratively Computed Eigenvalues",PETSC_DECIDE,PETSC_DECIDE,400,400,&ksp->eigviewer);
736: }
737: PetscViewerDrawGetDraw(ksp->eigviewer,0,&draw);
738: PetscDrawSPCreate(draw,1,&drawsp);
739: for (i=0; i<neig; i++) {
740: PetscDrawSPAddPoint(drawsp,r+i,c+i);
741: }
742: PetscDrawSPDraw(drawsp,PETSC_TRUE);
743: PetscDrawSPSave(drawsp);
744: PetscDrawSPDestroy(&drawsp);
745: }
746: if (flag3 && !rank) {
747: KSPPlotEigenContours_Private(ksp,neig,r,c);
748: }
749: PetscFree2(r,c);
750: }
751: }
753: PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_compute_singularvalues",NULL,NULL,&flag1);
754: if (flag1) {
755: PetscInt nits;
757: KSPGetIterationNumber(ksp,&nits);
758: if (!nits) {
759: PetscPrintf(comm,"Zero iterations in solver, cannot approximate any singular values\n");
760: } else {
761: PetscReal emax,emin;
763: KSPComputeExtremeSingularValues(ksp,&emax,&emin);
764: PetscPrintf(comm,"Iteratively computed extreme singular values: max %g min %g max/min %g\n",(double)emax,(double)emin,(double)(emax/emin));
765: }
766: }
768: PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_compute_eigenvalues_explicitly",NULL,NULL,&flag1);
769: PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_plot_eigenvalues_explicitly",NULL,NULL,&flag2);
770: if (flag1 || flag2) {
771: PetscInt n,i;
772: PetscReal *r,*c;
773: PetscMPIInt rank;
774: MPI_Comm_rank(comm,&rank);
775: VecGetSize(ksp->vec_sol,&n);
776: PetscMalloc2(n,&r,n,&c);
777: KSPComputeEigenvaluesExplicitly(ksp,n,r,c);
778: if (flag1) {
779: PetscPrintf(comm,"Explicitly computed eigenvalues\n");
780: for (i=0; i<n; i++) {
781: if (c[i] >= 0.0) {
782: PetscPrintf(comm,"%g + %gi\n",(double)r[i],(double)c[i]);
783: } else {
784: PetscPrintf(comm,"%g - %gi\n",(double)r[i],-(double)c[i]);
785: }
786: }
787: }
788: if (flag2 && !rank) {
789: PetscDraw draw;
790: PetscDrawSP drawsp;
792: if (!ksp->eigviewer) {
793: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Explicitly Computed Eigenvalues",PETSC_DECIDE,PETSC_DECIDE,400,400,&ksp->eigviewer);
794: }
795: PetscViewerDrawGetDraw(ksp->eigviewer,0,&draw);
796: PetscDrawSPCreate(draw,1,&drawsp);
797: PetscDrawSPReset(drawsp);
798: for (i=0; i<n; i++) {
799: PetscDrawSPAddPoint(drawsp,r+i,c+i);
800: }
801: PetscDrawSPDraw(drawsp,PETSC_TRUE);
802: PetscDrawSPSave(drawsp);
803: PetscDrawSPDestroy(&drawsp);
804: }
805: PetscFree2(r,c);
806: }
808: PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_view_mat_explicit",NULL,NULL,&flag2);
809: if (flag2) {
810: Mat A,B;
811: PCGetOperators(ksp->pc,&A,NULL);
812: MatComputeExplicitOperator(A,&B);
813: MatViewFromOptions(B,(PetscObject)ksp,"-ksp_view_mat_explicit");
814: MatDestroy(&B);
815: }
816: PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_view_preconditioned_operator_explicit",NULL,NULL,&flag2);
817: if (flag2) {
818: Mat B;
819: KSPComputeExplicitOperator(ksp,&B);
820: MatViewFromOptions(B,(PetscObject)ksp,"-ksp_view_preconditioned_operator_explicit");
821: MatDestroy(&B);
822: }
823: KSPViewFromOptions(ksp,NULL,"-ksp_view");
825: PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"-ksp_final_residual",NULL,NULL,&flg);
826: if (flg) {
827: Mat A;
828: Vec t;
829: PetscReal norm;
830: 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");
831: PCGetOperators(ksp->pc,&A,NULL);
832: VecDuplicate(ksp->vec_rhs,&t);
833: KSP_MatMult(ksp,A,ksp->vec_sol,t);
834: VecAYPX(t, -1.0, ksp->vec_rhs);
835: VecNorm(t,NORM_2,&norm);
836: VecDestroy(&t);
837: PetscPrintf(comm,"KSP final norm of residual %g\n",(double)norm);
838: }
839: VecViewFromOptions(ksp->vec_sol,(PetscObject)ksp,"-ksp_view_solution");
841: if (inXisinB) {
842: VecCopy(x,b);
843: VecDestroy(&x);
844: }
845: PetscObjectSAWsBlock((PetscObject)ksp);
846: if (ksp->errorifnotconverged && ksp->reason < 0) SETERRQ(comm,PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged");
847: return(0);
848: }
852: /*@
853: KSPSolveTranspose - Solves the transpose of a linear system.
855: Collective on KSP
857: Input Parameter:
858: + ksp - iterative context obtained from KSPCreate()
859: . b - right hand side vector
860: - x - solution vector
862: Notes: For complex numbers this solve the non-Hermitian transpose system.
864: Developer Notes: We need to implement a KSPSolveHermitianTranspose()
866: Level: developer
868: .keywords: KSP, solve, linear system
870: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPConvergedDefault(),
871: KSPSolve()
872: @*/
874: PetscErrorCode KSPSolveTranspose(KSP ksp,Vec b,Vec x)
875: {
877: PetscBool inXisinB=PETSC_FALSE;
883: if (x == b) {
884: VecDuplicate(b,&x);
885: inXisinB = PETSC_TRUE;
886: }
887: PetscObjectReference((PetscObject)b);
888: PetscObjectReference((PetscObject)x);
889: VecDestroy(&ksp->vec_rhs);
890: VecDestroy(&ksp->vec_sol);
892: ksp->vec_rhs = b;
893: ksp->vec_sol = x;
894: ksp->transpose_solve = PETSC_TRUE;
896: KSPSetUp(ksp);
897: if (ksp->guess_zero) { VecSet(ksp->vec_sol,0.0);}
898: (*ksp->ops->solve)(ksp);
899: if (!ksp->reason) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
900: KSPReasonViewFromOptions(ksp);
901: if (inXisinB) {
902: VecCopy(x,b);
903: VecDestroy(&x);
904: }
905: if (ksp->errorifnotconverged && ksp->reason < 0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged");
906: return(0);
907: }
911: /*@
912: KSPReset - Resets a KSP context to the kspsetupcalled = 0 state and removes any allocated Vecs and Mats
914: Collective on KSP
916: Input Parameter:
917: . ksp - iterative context obtained from KSPCreate()
919: Level: beginner
921: .keywords: KSP, destroy
923: .seealso: KSPCreate(), KSPSetUp(), KSPSolve()
924: @*/
925: PetscErrorCode KSPReset(KSP ksp)
926: {
931: if (!ksp) return(0);
932: if (ksp->ops->reset) {
933: (*ksp->ops->reset)(ksp);
934: }
935: if (ksp->pc) {PCReset(ksp->pc);}
936: KSPFischerGuessDestroy(&ksp->guess);
937: VecDestroyVecs(ksp->nwork,&ksp->work);
938: VecDestroy(&ksp->vec_rhs);
939: VecDestroy(&ksp->vec_sol);
940: VecDestroy(&ksp->diagonal);
941: VecDestroy(&ksp->truediagonal);
943: ksp->setupstage = KSP_SETUP_NEW;
944: return(0);
945: }
949: /*@
950: KSPDestroy - Destroys KSP context.
952: Collective on KSP
954: Input Parameter:
955: . ksp - iterative context obtained from KSPCreate()
957: Level: beginner
959: .keywords: KSP, destroy
961: .seealso: KSPCreate(), KSPSetUp(), KSPSolve()
962: @*/
963: PetscErrorCode KSPDestroy(KSP *ksp)
964: {
966: PC pc;
969: if (!*ksp) return(0);
971: if (--((PetscObject)(*ksp))->refct > 0) {*ksp = 0; return(0);}
973: PetscObjectSAWsViewOff((PetscObject)*ksp);
974: /*
975: Avoid a cascading call to PCReset(ksp->pc) from the following call:
976: PCReset() shouldn't be called from KSPDestroy() as it is unprotected by pc's
977: refcount (and may be shared, e.g., by other ksps).
978: */
979: pc = (*ksp)->pc;
980: (*ksp)->pc = NULL;
981: KSPReset((*ksp));
982: (*ksp)->pc = pc;
983: if ((*ksp)->ops->destroy) {(*(*ksp)->ops->destroy)(*ksp);}
985: DMDestroy(&(*ksp)->dm);
986: PCDestroy(&(*ksp)->pc);
987: PetscFree((*ksp)->res_hist_alloc);
988: if ((*ksp)->convergeddestroy) {
989: (*(*ksp)->convergeddestroy)((*ksp)->cnvP);
990: }
991: KSPMonitorCancel((*ksp));
992: PetscViewerDestroy(&(*ksp)->eigviewer);
993: PetscHeaderDestroy(ksp);
994: return(0);
995: }
999: /*@
1000: KSPSetPCSide - Sets the preconditioning side.
1002: Logically Collective on KSP
1004: Input Parameter:
1005: . ksp - iterative context obtained from KSPCreate()
1007: Output Parameter:
1008: . side - the preconditioning side, where side is one of
1009: .vb
1010: PC_LEFT - left preconditioning (default)
1011: PC_RIGHT - right preconditioning
1012: PC_SYMMETRIC - symmetric preconditioning
1013: .ve
1015: Options Database Keys:
1016: . -ksp_pc_side <right,left,symmetric>
1018: Notes:
1019: Left preconditioning is used by default for most Krylov methods except KSPFGMRES which only supports right preconditioning.
1021: For methods changing the side of the preconditioner changes the norm type that is used, see KSPSetNormType().
1023: Symmetric preconditioning is currently available only for the KSPQCG method. Note, however, that
1024: symmetric preconditioning can be emulated by using either right or left
1025: preconditioning and a pre or post processing step.
1027: Setting the PC side often affects the default norm type. See KSPSetNormType() for details.
1029: Level: intermediate
1031: .keywords: KSP, set, right, left, symmetric, side, preconditioner, flag
1033: .seealso: KSPGetPCSide(), KSPSetNormType(), KSPGetNormType()
1034: @*/
1035: PetscErrorCode KSPSetPCSide(KSP ksp,PCSide side)
1036: {
1040: ksp->pc_side = ksp->pc_side_set = side;
1041: return(0);
1042: }
1046: /*@
1047: KSPGetPCSide - Gets the preconditioning side.
1049: Not Collective
1051: Input Parameter:
1052: . ksp - iterative context obtained from KSPCreate()
1054: Output Parameter:
1055: . side - the preconditioning side, where side is one of
1056: .vb
1057: PC_LEFT - left preconditioning (default)
1058: PC_RIGHT - right preconditioning
1059: PC_SYMMETRIC - symmetric preconditioning
1060: .ve
1062: Level: intermediate
1064: .keywords: KSP, get, right, left, symmetric, side, preconditioner, flag
1066: .seealso: KSPSetPCSide()
1067: @*/
1068: PetscErrorCode KSPGetPCSide(KSP ksp,PCSide *side)
1069: {
1075: KSPSetUpNorms_Private(ksp,&ksp->normtype,&ksp->pc_side);
1076: *side = ksp->pc_side;
1077: return(0);
1078: }
1082: /*@
1083: KSPGetTolerances - Gets the relative, absolute, divergence, and maximum
1084: iteration tolerances used by the default KSP convergence tests.
1086: Not Collective
1088: Input Parameter:
1089: . ksp - the Krylov subspace context
1091: Output Parameters:
1092: + rtol - the relative convergence tolerance
1093: . abstol - the absolute convergence tolerance
1094: . dtol - the divergence tolerance
1095: - maxits - maximum number of iterations
1097: Notes:
1098: The user can specify NULL for any parameter that is not needed.
1100: Level: intermediate
1102: .keywords: KSP, get, tolerance, absolute, relative, divergence, convergence,
1103: maximum, iterations
1105: .seealso: KSPSetTolerances()
1106: @*/
1107: PetscErrorCode KSPGetTolerances(KSP ksp,PetscReal *rtol,PetscReal *abstol,PetscReal *dtol,PetscInt *maxits)
1108: {
1111: if (abstol) *abstol = ksp->abstol;
1112: if (rtol) *rtol = ksp->rtol;
1113: if (dtol) *dtol = ksp->divtol;
1114: if (maxits) *maxits = ksp->max_it;
1115: return(0);
1116: }
1120: /*@
1121: KSPSetTolerances - Sets the relative, absolute, divergence, and maximum
1122: iteration tolerances used by the default KSP convergence testers.
1124: Logically Collective on KSP
1126: Input Parameters:
1127: + ksp - the Krylov subspace context
1128: . rtol - the relative convergence tolerance, relative decrease in the (possibly preconditioned) residual norm
1129: . abstol - the absolute convergence tolerance absolute size of the (possibly preconditioned) residual norm
1130: . dtol - the divergence tolerance, amount (possibly preconditioned) residual norm can increase before KSPConvergedDefault() concludes that the method is diverging
1131: - maxits - maximum number of iterations to use
1133: Options Database Keys:
1134: + -ksp_atol <abstol> - Sets abstol
1135: . -ksp_rtol <rtol> - Sets rtol
1136: . -ksp_divtol <dtol> - Sets dtol
1137: - -ksp_max_it <maxits> - Sets maxits
1139: Notes:
1140: Use PETSC_DEFAULT to retain the default value of any of the tolerances.
1142: See KSPConvergedDefault() for details how these parameters are used in the default convergence test. See also KSPSetConvergenceTest()
1143: for setting user-defined stopping criteria.
1145: Level: intermediate
1147: .keywords: KSP, set, tolerance, absolute, relative, divergence,
1148: convergence, maximum, iterations
1150: .seealso: KSPGetTolerances(), KSPConvergedDefault(), KSPSetConvergenceTest()
1151: @*/
1152: PetscErrorCode KSPSetTolerances(KSP ksp,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt maxits)
1153: {
1161: if (rtol != PETSC_DEFAULT) {
1162: 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);
1163: ksp->rtol = rtol;
1164: }
1165: if (abstol != PETSC_DEFAULT) {
1166: if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
1167: ksp->abstol = abstol;
1168: }
1169: if (dtol != PETSC_DEFAULT) {
1170: if (dtol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Divergence tolerance %g must be larger than 1.0",(double)dtol);
1171: ksp->divtol = dtol;
1172: }
1173: if (maxits != PETSC_DEFAULT) {
1174: if (maxits < 0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxits);
1175: ksp->max_it = maxits;
1176: }
1177: return(0);
1178: }
1182: /*@
1183: KSPSetInitialGuessNonzero - Tells the iterative solver that the
1184: initial guess is nonzero; otherwise KSP assumes the initial guess
1185: is to be zero (and thus zeros it out before solving).
1187: Logically Collective on KSP
1189: Input Parameters:
1190: + ksp - iterative context obtained from KSPCreate()
1191: - flg - PETSC_TRUE indicates the guess is non-zero, PETSC_FALSE indicates the guess is zero
1193: Options database keys:
1194: . -ksp_initial_guess_nonzero : use nonzero initial guess; this takes an optional truth value (0/1/no/yes/true/false)
1196: Level: beginner
1198: Notes:
1199: If this is not called the X vector is zeroed in the call to KSPSolve().
1201: .keywords: KSP, set, initial guess, nonzero
1203: .seealso: KSPGetInitialGuessNonzero(), KSPSetInitialGuessKnoll(), KSPGetInitialGuessKnoll()
1204: @*/
1205: PetscErrorCode KSPSetInitialGuessNonzero(KSP ksp,PetscBool flg)
1206: {
1210: ksp->guess_zero = (PetscBool) !(int)flg;
1211: return(0);
1212: }
1216: /*@
1217: KSPGetInitialGuessNonzero - Determines whether the KSP solver is using
1218: a zero initial guess.
1220: Not Collective
1222: Input Parameter:
1223: . ksp - iterative context obtained from KSPCreate()
1225: Output Parameter:
1226: . flag - PETSC_TRUE if guess is nonzero, else PETSC_FALSE
1228: Level: intermediate
1230: .keywords: KSP, set, initial guess, nonzero
1232: .seealso: KSPSetInitialGuessNonzero(), KSPSetInitialGuessKnoll(), KSPGetInitialGuessKnoll()
1233: @*/
1234: PetscErrorCode KSPGetInitialGuessNonzero(KSP ksp,PetscBool *flag)
1235: {
1239: if (ksp->guess_zero) *flag = PETSC_FALSE;
1240: else *flag = PETSC_TRUE;
1241: return(0);
1242: }
1246: /*@
1247: KSPSetErrorIfNotConverged - Causes KSPSolve() to generate an error if the solver has not converged.
1249: Logically Collective on KSP
1251: Input Parameters:
1252: + ksp - iterative context obtained from KSPCreate()
1253: - flg - PETSC_TRUE indicates you want the error generated
1255: Options database keys:
1256: . -ksp_error_if_not_converged : this takes an optional truth value (0/1/no/yes/true/false)
1258: Level: intermediate
1260: Notes:
1261: Normally PETSc continues if a linear solver fails to converge, you can call KSPGetConvergedReason() after a KSPSolve()
1262: to determine if it has converged.
1264: .keywords: KSP, set, initial guess, nonzero
1266: .seealso: KSPGetInitialGuessNonzero(), KSPSetInitialGuessKnoll(), KSPGetInitialGuessKnoll(), KSPGetErrorIfNotConverged()
1267: @*/
1268: PetscErrorCode KSPSetErrorIfNotConverged(KSP ksp,PetscBool flg)
1269: {
1273: ksp->errorifnotconverged = flg;
1274: return(0);
1275: }
1279: /*@
1280: KSPGetErrorIfNotConverged - Will KSPSolve() generate an error if the solver does not converge?
1282: Not Collective
1284: Input Parameter:
1285: . ksp - iterative context obtained from KSPCreate()
1287: Output Parameter:
1288: . flag - PETSC_TRUE if it will generate an error, else PETSC_FALSE
1290: Level: intermediate
1292: .keywords: KSP, set, initial guess, nonzero
1294: .seealso: KSPSetInitialGuessNonzero(), KSPSetInitialGuessKnoll(), KSPGetInitialGuessKnoll(), KSPSetErrorIfNotConverged()
1295: @*/
1296: PetscErrorCode KSPGetErrorIfNotConverged(KSP ksp,PetscBool *flag)
1297: {
1301: *flag = ksp->errorifnotconverged;
1302: return(0);
1303: }
1307: /*@
1308: KSPSetInitialGuessKnoll - Tells the iterative solver to use PCApply(pc,b,..) to compute the initial guess (The Knoll trick)
1310: Logically Collective on KSP
1312: Input Parameters:
1313: + ksp - iterative context obtained from KSPCreate()
1314: - flg - PETSC_TRUE or PETSC_FALSE
1316: Level: advanced
1319: .keywords: KSP, set, initial guess, nonzero
1321: .seealso: KSPGetInitialGuessKnoll(), KSPSetInitialGuessNonzero(), KSPGetInitialGuessNonzero()
1322: @*/
1323: PetscErrorCode KSPSetInitialGuessKnoll(KSP ksp,PetscBool flg)
1324: {
1328: ksp->guess_knoll = flg;
1329: return(0);
1330: }
1334: /*@
1335: KSPGetInitialGuessKnoll - Determines whether the KSP solver is using the Knoll trick (using PCApply(pc,b,...) to compute
1336: the initial guess
1338: Not Collective
1340: Input Parameter:
1341: . ksp - iterative context obtained from KSPCreate()
1343: Output Parameter:
1344: . flag - PETSC_TRUE if using Knoll trick, else PETSC_FALSE
1346: Level: advanced
1348: .keywords: KSP, set, initial guess, nonzero
1350: .seealso: KSPSetInitialGuessKnoll(), KSPSetInitialGuessNonzero(), KSPGetInitialGuessNonzero()
1351: @*/
1352: PetscErrorCode KSPGetInitialGuessKnoll(KSP ksp,PetscBool *flag)
1353: {
1357: *flag = ksp->guess_knoll;
1358: return(0);
1359: }
1363: /*@
1364: KSPGetComputeSingularValues - Gets the flag indicating whether the extreme singular
1365: values will be calculated via a Lanczos or Arnoldi process as the linear
1366: system is solved.
1368: Not Collective
1370: Input Parameter:
1371: . ksp - iterative context obtained from KSPCreate()
1373: Output Parameter:
1374: . flg - PETSC_TRUE or PETSC_FALSE
1376: Options Database Key:
1377: . -ksp_monitor_singular_value - Activates KSPSetComputeSingularValues()
1379: Notes:
1380: Currently this option is not valid for all iterative methods.
1382: Many users may just want to use the monitoring routine
1383: KSPMonitorSingularValue() (which can be set with option -ksp_monitor_singular_value)
1384: to print the singular values at each iteration of the linear solve.
1386: Level: advanced
1388: .keywords: KSP, set, compute, singular values
1390: .seealso: KSPComputeExtremeSingularValues(), KSPMonitorSingularValue()
1391: @*/
1392: PetscErrorCode KSPGetComputeSingularValues(KSP ksp,PetscBool *flg)
1393: {
1397: *flg = ksp->calc_sings;
1398: return(0);
1399: }
1403: /*@
1404: KSPSetComputeSingularValues - Sets a flag so that the extreme singular
1405: values will be calculated via a Lanczos or Arnoldi process as the linear
1406: system is solved.
1408: Logically Collective on KSP
1410: Input Parameters:
1411: + ksp - iterative context obtained from KSPCreate()
1412: - flg - PETSC_TRUE or PETSC_FALSE
1414: Options Database Key:
1415: . -ksp_monitor_singular_value - Activates KSPSetComputeSingularValues()
1417: Notes:
1418: Currently this option is not valid for all iterative methods.
1420: Many users may just want to use the monitoring routine
1421: KSPMonitorSingularValue() (which can be set with option -ksp_monitor_singular_value)
1422: to print the singular values at each iteration of the linear solve.
1424: Level: advanced
1426: .keywords: KSP, set, compute, singular values
1428: .seealso: KSPComputeExtremeSingularValues(), KSPMonitorSingularValue()
1429: @*/
1430: PetscErrorCode KSPSetComputeSingularValues(KSP ksp,PetscBool flg)
1431: {
1435: ksp->calc_sings = flg;
1436: return(0);
1437: }
1441: /*@
1442: KSPGetComputeEigenvalues - Gets the flag indicating that the extreme eigenvalues
1443: values will be calculated via a Lanczos or Arnoldi process as the linear
1444: system is solved.
1446: Not Collective
1448: Input Parameter:
1449: . ksp - iterative context obtained from KSPCreate()
1451: Output Parameter:
1452: . flg - PETSC_TRUE or PETSC_FALSE
1454: Notes:
1455: Currently this option is not valid for all iterative methods.
1457: Level: advanced
1459: .keywords: KSP, set, compute, eigenvalues
1461: .seealso: KSPComputeEigenvalues(), KSPComputeEigenvaluesExplicitly()
1462: @*/
1463: PetscErrorCode KSPGetComputeEigenvalues(KSP ksp,PetscBool *flg)
1464: {
1468: *flg = ksp->calc_sings;
1469: return(0);
1470: }
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: }
1505: /*@
1506: KSPSetComputeRitz - Sets a flag so that the Ritz or harmonic Ritz pairs
1507: will be calculated via a Lanczos or Arnoldi process as the linear
1508: system is solved.
1510: Logically Collective on KSP
1512: Input Parameters:
1513: + ksp - iterative context obtained from KSPCreate()
1514: - flg - PETSC_TRUE or PETSC_FALSE
1516: Notes:
1517: Currently this option is only valid for the GMRES method.
1519: Level: advanced
1521: .keywords: KSP, set, compute, ritz
1523: .seealso: KSPComputeRitz()
1524: @*/
1525: PetscErrorCode KSPSetComputeRitz(KSP ksp, PetscBool flg)
1526: {
1530: ksp->calc_ritz = flg;
1531: return(0);
1532: }
1536: /*@
1537: KSPGetRhs - Gets the right-hand-side vector for the linear system to
1538: be solved.
1540: Not Collective
1542: Input Parameter:
1543: . ksp - iterative context obtained from KSPCreate()
1545: Output Parameter:
1546: . r - right-hand-side vector
1548: Level: developer
1550: .keywords: KSP, get, right-hand-side, rhs
1552: .seealso: KSPGetSolution(), KSPSolve()
1553: @*/
1554: PetscErrorCode KSPGetRhs(KSP ksp,Vec *r)
1555: {
1559: *r = ksp->vec_rhs;
1560: return(0);
1561: }
1565: /*@
1566: KSPGetSolution - Gets the location of the solution for the
1567: linear system to be solved. Note that this may not be where the solution
1568: is stored during the iterative process; see KSPBuildSolution().
1570: Not Collective
1572: Input Parameters:
1573: . ksp - iterative context obtained from KSPCreate()
1575: Output Parameters:
1576: . v - solution vector
1578: Level: developer
1580: .keywords: KSP, get, solution
1582: .seealso: KSPGetRhs(), KSPBuildSolution(), KSPSolve()
1583: @*/
1584: PetscErrorCode KSPGetSolution(KSP ksp,Vec *v)
1585: {
1589: *v = ksp->vec_sol;
1590: return(0);
1591: }
1595: /*@
1596: KSPSetPC - Sets the preconditioner to be used to calculate the
1597: application of the preconditioner on a vector.
1599: Collective on KSP
1601: Input Parameters:
1602: + ksp - iterative context obtained from KSPCreate()
1603: - pc - the preconditioner object
1605: Notes:
1606: Use KSPGetPC() to retrieve the preconditioner context (for example,
1607: to free it at the end of the computations).
1609: Level: developer
1611: .keywords: KSP, set, precondition, Binv
1613: .seealso: KSPGetPC()
1614: @*/
1615: PetscErrorCode KSPSetPC(KSP ksp,PC pc)
1616: {
1623: PetscObjectReference((PetscObject)pc);
1624: PCDestroy(&ksp->pc);
1625: ksp->pc = pc;
1626: PetscLogObjectParent((PetscObject)ksp,(PetscObject)ksp->pc);
1627: return(0);
1628: }
1632: /*@
1633: KSPGetPC - Returns a pointer to the preconditioner context
1634: set with KSPSetPC().
1636: Not Collective
1638: Input Parameters:
1639: . ksp - iterative context obtained from KSPCreate()
1641: Output Parameter:
1642: . pc - preconditioner context
1644: Level: developer
1646: .keywords: KSP, get, preconditioner, Binv
1648: .seealso: KSPSetPC()
1649: @*/
1650: PetscErrorCode KSPGetPC(KSP ksp,PC *pc)
1651: {
1657: if (!ksp->pc) {
1658: PCCreate(PetscObjectComm((PetscObject)ksp),&ksp->pc);
1659: PetscObjectIncrementTabLevel((PetscObject)ksp->pc,(PetscObject)ksp,0);
1660: PetscLogObjectParent((PetscObject)ksp,(PetscObject)ksp->pc);
1661: }
1662: *pc = ksp->pc;
1663: return(0);
1664: }
1668: /*@
1669: KSPMonitor - runs the user provided monitor routines, if they exist
1671: Collective on KSP
1673: Input Parameters:
1674: + ksp - iterative context obtained from KSPCreate()
1675: . it - iteration number
1676: - rnorm - relative norm of the residual
1678: Notes:
1679: This routine is called by the KSP implementations.
1680: It does not typically need to be called by the user.
1682: Level: developer
1684: .seealso: KSPMonitorSet()
1685: @*/
1686: PetscErrorCode KSPMonitor(KSP ksp,PetscInt it,PetscReal rnorm)
1687: {
1688: PetscInt i, n = ksp->numbermonitors;
1692: for (i=0; i<n; i++) {
1693: (*ksp->monitor[i])(ksp,it,rnorm,ksp->monitorcontext[i]);
1694: }
1695: return(0);
1696: }
1700: /*@C
1701: KSPMonitorSet - Sets an ADDITIONAL function to be called at every iteration to monitor
1702: the residual/error etc.
1704: Logically Collective on KSP
1706: Input Parameters:
1707: + ksp - iterative context obtained from KSPCreate()
1708: . monitor - pointer to function (if this is NULL, it turns off monitoring
1709: . mctx - [optional] context for private data for the
1710: monitor routine (use NULL if no context is desired)
1711: - monitordestroy - [optional] routine that frees monitor context
1712: (may be NULL)
1714: Calling Sequence of monitor:
1715: $ monitor (KSP ksp, int it, PetscReal rnorm, void *mctx)
1717: + ksp - iterative context obtained from KSPCreate()
1718: . it - iteration number
1719: . rnorm - (estimated) 2-norm of (preconditioned) residual
1720: - mctx - optional monitoring context, as set by KSPMonitorSet()
1722: Options Database Keys:
1723: + -ksp_monitor - sets KSPMonitorDefault()
1724: . -ksp_monitor_true_residual - sets KSPMonitorTrueResidualNorm()
1725: . -ksp_monitor_max - sets KSPMonitorTrueResidualMaxNorm()
1726: . -ksp_monitor_lg_residualnorm - sets line graph monitor,
1727: uses KSPMonitorLGResidualNormCreate()
1728: . -ksp_monitor_lg_true_residualnorm - sets line graph monitor,
1729: uses KSPMonitorLGResidualNormCreate()
1730: . -ksp_monitor_singular_value - sets KSPMonitorSingularValue()
1731: - -ksp_monitor_cancel - cancels all monitors that have
1732: been hardwired into a code by
1733: calls to KSPMonitorSet(), but
1734: does not cancel those set via
1735: the options database.
1737: Notes:
1738: The default is to do nothing. To print the residual, or preconditioned
1739: residual if KSPSetNormType(ksp,KSP_NORM_PRECONDITIONED) was called, use
1740: KSPMonitorDefault() as the monitoring routine, with a ASCII viewer as the
1741: context.
1743: Several different monitoring routines may be set by calling
1744: KSPMonitorSet() multiple times; all will be called in the
1745: order in which they were set.
1747: Fortran notes: Only a single monitor function can be set for each KSP object
1749: Level: beginner
1751: .keywords: KSP, set, monitor
1753: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(), KSPMonitorCancel()
1754: @*/
1755: PetscErrorCode KSPMonitorSet(KSP ksp,PetscErrorCode (*monitor)(KSP,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
1756: {
1757: PetscInt i;
1762: if (ksp->numbermonitors >= MAXKSPMONITORS) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Too many KSP monitors set");
1763: for (i=0; i<ksp->numbermonitors;i++) {
1764: if (monitor == ksp->monitor[i] && monitordestroy == ksp->monitordestroy[i] && mctx == ksp->monitorcontext[i]) {
1765: if (monitordestroy) {
1766: (*monitordestroy)(&mctx);
1767: }
1768: return(0);
1769: }
1770: }
1771: ksp->monitor[ksp->numbermonitors] = monitor;
1772: ksp->monitordestroy[ksp->numbermonitors] = monitordestroy;
1773: ksp->monitorcontext[ksp->numbermonitors++] = (void*)mctx;
1774: return(0);
1775: }
1779: /*@
1780: KSPMonitorCancel - Clears all monitors for a KSP object.
1782: Logically Collective on KSP
1784: Input Parameters:
1785: . ksp - iterative context obtained from KSPCreate()
1787: Options Database Key:
1788: . -ksp_monitor_cancel - Cancels all monitors that have
1789: been hardwired into a code by calls to KSPMonitorSet(),
1790: but does not cancel those set via the options database.
1792: Level: intermediate
1794: .keywords: KSP, set, monitor
1796: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(), KSPMonitorSet()
1797: @*/
1798: PetscErrorCode KSPMonitorCancel(KSP ksp)
1799: {
1801: PetscInt i;
1805: for (i=0; i<ksp->numbermonitors; i++) {
1806: if (ksp->monitordestroy[i]) {
1807: (*ksp->monitordestroy[i])(&ksp->monitorcontext[i]);
1808: }
1809: }
1810: ksp->numbermonitors = 0;
1811: return(0);
1812: }
1816: /*@C
1817: KSPGetMonitorContext - Gets the monitoring context, as set by
1818: KSPMonitorSet() for the FIRST monitor only.
1820: Not Collective
1822: Input Parameter:
1823: . ksp - iterative context obtained from KSPCreate()
1825: Output Parameter:
1826: . ctx - monitoring context
1828: Level: intermediate
1830: .keywords: KSP, get, monitor, context
1832: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate()
1833: @*/
1834: PetscErrorCode KSPGetMonitorContext(KSP ksp,void **ctx)
1835: {
1838: *ctx = (ksp->monitorcontext[0]);
1839: return(0);
1840: }
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: }
1896: /*@C
1897: KSPGetResidualHistory - Gets the array used to hold the residual history
1898: and the number of residuals it contains.
1900: Not Collective
1902: Input Parameter:
1903: . ksp - iterative context obtained from KSPCreate()
1905: Output Parameters:
1906: + a - pointer to array to hold history (or NULL)
1907: - na - number of used entries in a (or NULL)
1909: Level: advanced
1911: Notes:
1912: Can only be called after a KSPSetResidualHistory() otherwise a and na are set to zero
1914: The Fortran version of this routine has a calling sequence
1915: $ call KSPGetResidualHistory(KSP ksp, integer na, integer ierr)
1916: note that you have passed a Fortran array into KSPSetResidualHistory() and you need
1917: to access the residual values from this Fortran array you provided. Only the na (number of
1918: residual norms currently held) is set.
1920: .keywords: KSP, get, residual, history, norm
1922: .seealso: KSPGetResidualHistory()
1924: @*/
1925: PetscErrorCode KSPGetResidualHistory(KSP ksp,PetscReal *a[],PetscInt *na)
1926: {
1929: if (a) *a = ksp->res_hist;
1930: if (na) *na = ksp->res_hist_len;
1931: return(0);
1932: }
1936: /*@C
1937: KSPSetConvergenceTest - Sets the function to be used to determine
1938: convergence.
1940: Logically Collective on KSP
1942: Input Parameters:
1943: + ksp - iterative context obtained from KSPCreate()
1944: . converge - pointer to int function
1945: . cctx - context for private data for the convergence routine (may be null)
1946: - destroy - a routine for destroying the context (may be null)
1948: Calling sequence of converge:
1949: $ converge (KSP ksp, int it, PetscReal rnorm, KSPConvergedReason *reason,void *mctx)
1951: + ksp - iterative context obtained from KSPCreate()
1952: . it - iteration number
1953: . rnorm - (estimated) 2-norm of (preconditioned) residual
1954: . reason - the reason why it has converged or diverged
1955: - cctx - optional convergence context, as set by KSPSetConvergenceTest()
1958: Notes:
1959: Must be called after the KSP type has been set so put this after
1960: a call to KSPSetType(), or KSPSetFromOptions().
1962: The default convergence test, KSPConvergedDefault(), aborts if the
1963: residual grows to more than 10000 times the initial residual.
1965: The default is a combination of relative and absolute tolerances.
1966: The residual value that is tested may be an approximation; routines
1967: that need exact values should compute them.
1969: In the default PETSc convergence test, the precise values of reason
1970: are macros such as KSP_CONVERGED_RTOL, which are defined in petscksp.h.
1972: Level: advanced
1974: .keywords: KSP, set, convergence, test, context
1976: .seealso: KSPConvergedDefault(), KSPGetConvergenceContext(), KSPSetTolerances()
1977: @*/
1978: PetscErrorCode KSPSetConvergenceTest(KSP ksp,PetscErrorCode (*converge)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
1979: {
1984: if (ksp->convergeddestroy) {
1985: (*ksp->convergeddestroy)(ksp->cnvP);
1986: }
1987: ksp->converged = converge;
1988: ksp->convergeddestroy = destroy;
1989: ksp->cnvP = (void*)cctx;
1990: return(0);
1991: }
1995: /*@C
1996: KSPGetConvergenceContext - Gets the convergence context set with
1997: KSPSetConvergenceTest().
1999: Not Collective
2001: Input Parameter:
2002: . ksp - iterative context obtained from KSPCreate()
2004: Output Parameter:
2005: . ctx - monitoring context
2007: Level: advanced
2009: .keywords: KSP, get, convergence, test, context
2011: .seealso: KSPConvergedDefault(), KSPSetConvergenceTest()
2012: @*/
2013: PetscErrorCode KSPGetConvergenceContext(KSP ksp,void **ctx)
2014: {
2017: *ctx = ksp->cnvP;
2018: return(0);
2019: }
2023: /*@C
2024: KSPBuildSolution - Builds the approximate solution in a vector provided.
2025: This routine is NOT commonly needed (see KSPSolve()).
2027: Collective on KSP
2029: Input Parameter:
2030: . ctx - iterative context obtained from KSPCreate()
2032: Output Parameter:
2033: Provide exactly one of
2034: + v - location to stash solution.
2035: - V - the solution is returned in this location. This vector is created
2036: internally. This vector should NOT be destroyed by the user with
2037: VecDestroy().
2039: Notes:
2040: This routine can be used in one of two ways
2041: .vb
2042: KSPBuildSolution(ksp,NULL,&V);
2043: or
2044: KSPBuildSolution(ksp,v,NULL); or KSPBuildSolution(ksp,v,&v);
2045: .ve
2046: In the first case an internal vector is allocated to store the solution
2047: (the user cannot destroy this vector). In the second case the solution
2048: is generated in the vector that the user provides. Note that for certain
2049: methods, such as KSPCG, the second case requires a copy of the solution,
2050: while in the first case the call is essentially free since it simply
2051: returns the vector where the solution already is stored. For some methods
2052: like GMRES this is a reasonably expensive operation and should only be
2053: used in truly needed.
2055: Level: advanced
2057: .keywords: KSP, build, solution
2059: .seealso: KSPGetSolution(), KSPBuildResidual()
2060: @*/
2061: PetscErrorCode KSPBuildSolution(KSP ksp,Vec v,Vec *V)
2062: {
2067: if (!V && !v) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONG,"Must provide either v or V");
2068: if (!V) V = &v;
2069: (*ksp->ops->buildsolution)(ksp,v,V);
2070: return(0);
2071: }
2075: /*@C
2076: KSPBuildResidual - Builds the residual in a vector provided.
2078: Collective on KSP
2080: Input Parameter:
2081: . ksp - iterative context obtained from KSPCreate()
2083: Output Parameters:
2084: + v - optional location to stash residual. If v is not provided,
2085: then a location is generated.
2086: . t - work vector. If not provided then one is generated.
2087: - V - the residual
2089: Notes:
2090: Regardless of whether or not v is provided, the residual is
2091: returned in V.
2093: Level: advanced
2095: .keywords: KSP, build, residual
2097: .seealso: KSPBuildSolution()
2098: @*/
2099: PetscErrorCode KSPBuildResidual(KSP ksp,Vec t,Vec v,Vec *V)
2100: {
2102: PetscBool flag = PETSC_FALSE;
2103: Vec w = v,tt = t;
2107: if (!w) {
2108: VecDuplicate(ksp->vec_rhs,&w);
2109: PetscLogObjectParent((PetscObject)ksp,(PetscObject)w);
2110: }
2111: if (!tt) {
2112: VecDuplicate(ksp->vec_sol,&tt); flag = PETSC_TRUE;
2113: PetscLogObjectParent((PetscObject)ksp,(PetscObject)tt);
2114: }
2115: (*ksp->ops->buildresidual)(ksp,tt,w,V);
2116: if (flag) {VecDestroy(&tt);}
2117: return(0);
2118: }
2122: /*@
2123: KSPSetDiagonalScale - Tells KSP to symmetrically diagonally scale the system
2124: before solving. This actually CHANGES the matrix (and right hand side).
2126: Logically Collective on KSP
2128: Input Parameter:
2129: + ksp - the KSP context
2130: - scale - PETSC_TRUE or PETSC_FALSE
2132: Options Database Key:
2133: + -ksp_diagonal_scale -
2134: - -ksp_diagonal_scale_fix - scale the matrix back AFTER the solve
2137: Notes: Scales the matrix by D^(-1/2) A D^(-1/2) [D^(1/2) x ] = D^(-1/2) b
2138: where D_{ii} is 1/abs(A_{ii}) unless A_{ii} is zero and then it is 1.
2140: BE CAREFUL with this routine: it actually scales the matrix and right
2141: hand side that define the system. After the system is solved the matrix
2142: and right hand side remain scaled unless you use KSPSetDiagonalScaleFix()
2144: This should NOT be used within the SNES solves if you are using a line
2145: search.
2147: If you use this with the PCType Eisenstat preconditioner than you can
2148: use the PCEisenstatSetNoDiagonalScaling() option, or -pc_eisenstat_no_diagonal_scaling
2149: to save some unneeded, redundant flops.
2151: Level: intermediate
2153: .keywords: KSP, set, options, prefix, database
2155: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScaleFix()
2156: @*/
2157: PetscErrorCode KSPSetDiagonalScale(KSP ksp,PetscBool scale)
2158: {
2162: ksp->dscale = scale;
2163: return(0);
2164: }
2168: /*@
2169: KSPGetDiagonalScale - Checks if KSP solver scales the matrix and
2170: right hand side
2172: Not Collective
2174: Input Parameter:
2175: . ksp - the KSP context
2177: Output Parameter:
2178: . scale - PETSC_TRUE or PETSC_FALSE
2180: Notes:
2181: BE CAREFUL with this routine: it actually scales the matrix and right
2182: hand side that define the system. After the system is solved the matrix
2183: and right hand side remain scaled unless you use KSPSetDiagonalScaleFix()
2185: Level: intermediate
2187: .keywords: KSP, set, options, prefix, database
2189: .seealso: KSPSetDiagonalScale(), KSPSetDiagonalScaleFix()
2190: @*/
2191: PetscErrorCode KSPGetDiagonalScale(KSP ksp,PetscBool *scale)
2192: {
2196: *scale = ksp->dscale;
2197: return(0);
2198: }
2202: /*@
2203: KSPSetDiagonalScaleFix - Tells KSP to diagonally scale the system
2204: back after solving.
2206: Logically Collective on KSP
2208: Input Parameter:
2209: + ksp - the KSP context
2210: - fix - PETSC_TRUE to scale back after the system solve, PETSC_FALSE to not
2211: rescale (default)
2213: Notes:
2214: Must be called after KSPSetDiagonalScale()
2216: Using this will slow things down, because it rescales the matrix before and
2217: after each linear solve. This is intended mainly for testing to allow one
2218: to easily get back the original system to make sure the solution computed is
2219: accurate enough.
2221: Level: intermediate
2223: .keywords: KSP, set, options, prefix, database
2225: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScale(), KSPGetDiagonalScaleFix()
2226: @*/
2227: PetscErrorCode KSPSetDiagonalScaleFix(KSP ksp,PetscBool fix)
2228: {
2232: ksp->dscalefix = fix;
2233: return(0);
2234: }
2238: /*@
2239: KSPGetDiagonalScaleFix - Determines if KSP diagonally scales the system
2240: back after solving.
2242: Not Collective
2244: Input Parameter:
2245: . ksp - the KSP context
2247: Output Parameter:
2248: . fix - PETSC_TRUE to scale back after the system solve, PETSC_FALSE to not
2249: rescale (default)
2251: Notes:
2252: Must be called after KSPSetDiagonalScale()
2254: If PETSC_TRUE will slow things down, because it rescales the matrix before and
2255: after each linear solve. This is intended mainly for testing to allow one
2256: to easily get back the original system to make sure the solution computed is
2257: accurate enough.
2259: Level: intermediate
2261: .keywords: KSP, set, options, prefix, database
2263: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScale(), KSPSetDiagonalScaleFix()
2264: @*/
2265: PetscErrorCode KSPGetDiagonalScaleFix(KSP ksp,PetscBool *fix)
2266: {
2270: *fix = ksp->dscalefix;
2271: return(0);
2272: }
2276: /*@C
2277: KSPSetComputeOperators - set routine to compute the linear operators
2279: Logically Collective
2281: Input Arguments:
2282: + ksp - the KSP context
2283: . func - function to compute the operators
2284: - ctx - optional context
2286: Calling sequence of func:
2287: $ func(KSP ksp,Mat A,Mat B,void *ctx)
2289: + ksp - the KSP context
2290: . A - the linear operator
2291: . B - preconditioning matrix
2292: - ctx - optional user-provided context
2294: 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
2295: unless either KSPSetComputeOperators() or KSPSetOperators() is called before that KSPSolve() is called.
2297: To reuse the same preconditioner for the next KSPSolve() and not compute a new one based on the most recently computed matrix call KSPSetReusePreconditioner()
2299: Level: beginner
2301: .seealso: KSPSetOperators(), KSPSetComputeRHS(), DMKSPSetComputeOperators(), KSPSetComputeInitialGuess()
2302: @*/
2303: PetscErrorCode KSPSetComputeOperators(KSP ksp,PetscErrorCode (*func)(KSP,Mat,Mat,void*),void *ctx)
2304: {
2306: DM dm;
2310: KSPGetDM(ksp,&dm);
2311: DMKSPSetComputeOperators(dm,func,ctx);
2312: if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX;
2313: return(0);
2314: }
2318: /*@C
2319: KSPSetComputeRHS - set routine to compute the right hand side of the linear system
2321: Logically Collective
2323: Input Arguments:
2324: + ksp - the KSP context
2325: . func - function to compute the right hand side
2326: - ctx - optional context
2328: Calling sequence of func:
2329: $ func(KSP ksp,Vec b,void *ctx)
2331: + ksp - the KSP context
2332: . b - right hand side of linear system
2333: - ctx - optional user-provided context
2335: Notes: The routine you provide will be called EACH you call KSPSolve() to prepare the new right hand side for that solve
2337: Level: beginner
2339: .seealso: KSPSolve(), DMKSPSetComputeRHS(), KSPSetComputeOperators()
2340: @*/
2341: PetscErrorCode KSPSetComputeRHS(KSP ksp,PetscErrorCode (*func)(KSP,Vec,void*),void *ctx)
2342: {
2344: DM dm;
2348: KSPGetDM(ksp,&dm);
2349: DMKSPSetComputeRHS(dm,func,ctx);
2350: return(0);
2351: }
2355: /*@C
2356: KSPSetComputeInitialGuess - set routine to compute the initial guess of the linear system
2358: Logically Collective
2360: Input Arguments:
2361: + ksp - the KSP context
2362: . func - function to compute the initial guess
2363: - ctx - optional context
2365: Calling sequence of func:
2366: $ func(KSP ksp,Vec x,void *ctx)
2368: + ksp - the KSP context
2369: . x - solution vector
2370: - ctx - optional user-provided context
2372: Level: beginner
2374: .seealso: KSPSolve(), KSPSetComputeRHS(), KSPSetComputeOperators(), DMKSPSetComputeInitialGuess()
2375: @*/
2376: PetscErrorCode KSPSetComputeInitialGuess(KSP ksp,PetscErrorCode (*func)(KSP,Vec,void*),void *ctx)
2377: {
2379: DM dm;
2383: KSPGetDM(ksp,&dm);
2384: DMKSPSetComputeInitialGuess(dm,func,ctx);
2385: return(0);
2386: }