Actual source code: itfunc.c
1: /*$Id: itfunc.c,v 1.159 2001/08/07 03:03:45 balay Exp $*/
2: /*
3: Interface KSP routines that the user calls.
4: */
6: #include src/sles/ksp/kspimpl.h
8: /*@
9: KSPComputeExtremeSingularValues - Computes the extreme singular values
10: for the preconditioned operator. Called after or during KSPSolve()
11: (SLESSolve()).
13: Not Collective
15: Input Parameter:
16: . ksp - iterative context obtained from KSPCreate()
18: Output Parameters:
19: . emin, emax - extreme singular values
21: Notes:
22: One must call KSPSetComputeSingularValues() before calling KSPSetUp()
23: (or use the option -ksp_compute_eigenvalues) in order for this routine to work correctly.
25: Many users may just want to use the monitoring routine
26: KSPSingularValueMonitor() (which can be set with option -ksp_singmonitor)
27: to print the extreme singular values at each iteration of the linear solve.
29: Level: advanced
31: .keywords: KSP, compute, extreme, singular, values
33: .seealso: KSPSetComputeSingularValues(), KSPSingularValueMonitor(), KSPComputeEigenvalues()
34: @*/
35: int KSPComputeExtremeSingularValues(KSP ksp,PetscReal *emax,PetscReal *emin)
36: {
43: if (!ksp->calc_sings) {
44: SETERRQ(4,"Singular values not requested before KSPSetUp()");
45: }
47: if (ksp->ops->computeextremesingularvalues) {
48: (*ksp->ops->computeextremesingularvalues)(ksp,emax,emin);
49: }
50: return(0);
51: }
53: /*@
54: KSPComputeEigenvalues - Computes the extreme eigenvalues for the
55: preconditioned operator. Called after or during KSPSolve() (SLESSolve()).
57: Not Collective
59: Input Parameter:
60: + ksp - iterative context obtained from KSPCreate()
61: - n - size of arrays r and c. The number of eigenvalues computed (neig) will, in
62: general, be less than this.
64: Output Parameters:
65: + r - real part of computed eigenvalues
66: . c - complex part of computed eigenvalues
67: - neig - number of eigenvalues computed (will be less than or equal to n)
69: Options Database Keys:
70: + -ksp_compute_eigenvalues - Prints eigenvalues to stdout
71: - -ksp_plot_eigenvalues - Plots eigenvalues in an x-window display
73: Notes:
74: The number of eigenvalues estimated depends on the size of the Krylov space
75: generated during the KSPSolve() (that is the SLESSolve); for example, with
76: CG it corresponds to the number of CG iterations, for GMRES it is the number
77: of GMRES iterations SINCE the last restart. Any extra space in r[] and c[]
78: will be ignored.
80: KSPComputeEigenvalues() does not usually provide accurate estimates; it is
81: intended only for assistance in understanding the convergence of iterative
82: methods, not for eigenanalysis.
84: One must call KSPSetComputeEigenvalues() before calling KSPSetUp()
85: in order for this routine to work correctly.
87: Many users may just want to use the monitoring routine
88: KSPSingularValueMonitor() (which can be set with option -ksp_singmonitor)
89: to print the singular values at each iteration of the linear solve.
91: Level: advanced
93: .keywords: KSP, compute, extreme, singular, values
95: .seealso: KSPSetComputeSingularValues(), KSPSingularValueMonitor(), KSPComputeExtremeSingularValues()
96: @*/
97: int KSPComputeEigenvalues(KSP ksp,int n,PetscReal *r,PetscReal *c,int *neig)
98: {
105: if (!ksp->calc_sings) {
106: SETERRQ(4,"Eigenvalues not requested before KSPSetUp()");
107: }
109: if (ksp->ops->computeeigenvalues) {
110: (*ksp->ops->computeeigenvalues)(ksp,n,r,c,neig);
111: }
112: return(0);
113: }
115: /*@
116: KSPSetUp - Sets up the internal data structures for the
117: later use of an iterative solver.
119: Collective on KSP
121: Input Parameter:
122: . ksp - iterative context obtained from KSPCreate()
124: Level: developer
126: .keywords: KSP, setup
128: .seealso: KSPCreate(), KSPSolve(), KSPDestroy()
129: @*/
130: int KSPSetUp(KSP ksp)
131: {
137: /* reset the convergence flag from the previous solves */
138: ksp->reason = KSP_CONVERGED_ITERATING;
140: if (!ksp->type_name){
141: KSPSetType(ksp,KSPGMRES);
142: }
144: if (ksp->setupcalled) return(0);
145: ksp->setupcalled = 1;
146: (*ksp->ops->setup)(ksp);
147: return(0);
148: }
151: /*@
152: KSPSolve - Solves linear system; usually not called directly, rather
153: it is called by a call to SLESSolve().
155: Collective on KSP
157: Input Parameter:
158: . ksp - iterative context obtained from KSPCreate()
160: Output Parameter:
161: . its - number of iterations required
163: Options Database:
164: + -ksp_compute_eigenvalues - compute preconditioned operators eigenvalues
165: . -ksp_plot_eigenvalues - plot the computed eigenvalues in an X-window
166: . -ksp_compute_eigenvalues_explicitly - compute the eigenvalues by forming the
167: dense operator and useing LAPACK
168: - -ksp_plot_eigenvalues_explicitly - plot the explicitly computing eigenvalues
170: Notes:
171: On return, the parameter "its" contains either the iteration
172: number at which convergence was successfully reached or failure was detected.
174: Call KSPGetConvergedReason() to determine if the solver converged or failed and
175: why.
176:
177: If using a direct method (e.g., via the KSP solver
178: KSPPREONLY and a preconditioner such as PCLU/PCILU),
179: then its=1. See KSPSetTolerances() and KSPDefaultConverged()
180: for more details.
182: Understanding Convergence:
183: The routines KSPSetMonitor(), KSPComputeEigenvalues(), and
184: KSPComputeEigenvaluesExplicitly() provide information on additional
185: options to monitor convergence and print eigenvalue information.
187: Level: developer
189: .keywords: KSP, solve, linear system
191: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPDefaultConverged(),
192: SLESSolve(), KSPSolveTranspose(), SLESGetKSP()
193: @*/
194: int KSPSolve(KSP ksp,int *its)
195: {
196: int ierr,rank,nits;
197: PetscTruth flag1,flag2;
198: PetscScalar zero = 0.0;
203: if (!ksp->setupcalled){ KSPSetUp(ksp);}
204: if (ksp->guess_zero) { VecSet(&zero,ksp->vec_sol);}
205: /* reset the residual history list if requested */
206: if (ksp->res_hist_reset) ksp->res_hist_len = 0;
208: ksp->transpose_solve = PETSC_FALSE;
209: (*ksp->ops->solve)(ksp,&nits);
210: if (!ksp->reason) {
211: SETERRQ(1,"Internal error, solver returned without setting converged reason");
212: }
213: if (its) *its = nits;
215: MPI_Comm_rank(ksp->comm,&rank);
217: PetscOptionsHasName(ksp->prefix,"-ksp_compute_eigenvalues",&flag1);
218: PetscOptionsHasName(ksp->prefix,"-ksp_plot_eigenvalues",&flag2);
219: if (flag1 || flag2) {
220: int n = nits + 2,i,neig;
221: PetscReal *r,*c;
223: if (!n) {
224: PetscPrintf(ksp->comm,"Zero iterations in solver, cannot approximate any eigenvaluesn");
225: } else {
226: PetscMalloc(2*n*sizeof(PetscReal),&r);
227: c = r + n;
228: KSPComputeEigenvalues(ksp,n,r,c,&neig);
229: if (flag1) {
230: PetscPrintf(ksp->comm,"Iteratively computed eigenvaluesn");
231: for (i=0; i<neig; i++) {
232: if (c[i] >= 0.0) {PetscPrintf(ksp->comm,"%g + %gin",r[i],c[i]);}
233: else {PetscPrintf(ksp->comm,"%g - %gin",r[i],-c[i]);}
234: }
235: }
236: if (flag2 && !rank) {
237: PetscViewer viewer;
238: PetscDraw draw;
239: PetscDrawSP drawsp;
241: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Iteratively Computed Eigenvalues",
242: PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer);
243: PetscViewerDrawGetDraw(viewer,0,&draw);
244: PetscDrawSPCreate(draw,1,&drawsp);
245: for (i=0; i<neig; i++) {
246: PetscDrawSPAddPoint(drawsp,r+i,c+i);
247: }
248: PetscDrawSPDraw(drawsp);
249: PetscDrawSPDestroy(drawsp);
250: PetscViewerDestroy(viewer);
251: }
252: PetscFree(r);
253: }
254: }
256: PetscOptionsHasName(ksp->prefix,"-ksp_compute_eigenvalues_explicitly",&flag1);
257: PetscOptionsHasName(ksp->prefix,"-ksp_plot_eigenvalues_explicitly",&flag2);
258: if (flag1 || flag2) {
259: int n,i;
260: PetscReal *r,*c;
261: VecGetSize(ksp->vec_sol,&n);
262: PetscMalloc(2*n*sizeof(PetscReal),&r);
263: c = r + n;
264: KSPComputeEigenvaluesExplicitly(ksp,n,r,c);
265: if (flag1) {
266: PetscPrintf(ksp->comm,"Explicitly computed eigenvaluesn");
267: for (i=0; i<n; i++) {
268: if (c[i] >= 0.0) {PetscPrintf(ksp->comm,"%g + %gin",r[i],c[i]);}
269: else {PetscPrintf(ksp->comm,"%g - %gin",r[i],-c[i]);}
270: }
271: }
272: if (flag2 && !rank) {
273: PetscViewer viewer;
274: PetscDraw draw;
275: PetscDrawSP drawsp;
277: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Explicitly Computed Eigenvalues",0,320,300,300,&viewer);
278: PetscViewerDrawGetDraw(viewer,0,&draw);
279: PetscDrawSPCreate(draw,1,&drawsp);
280: for (i=0; i<n; i++) {
281: PetscDrawSPAddPoint(drawsp,r+i,c+i);
282: }
283: PetscDrawSPDraw(drawsp);
284: PetscDrawSPDestroy(drawsp);
285: PetscViewerDestroy(viewer);
286: }
287: PetscFree(r);
288: }
290: PetscOptionsHasName(ksp->prefix,"-ksp_view_operator",&flag2);
291: if (flag2) {
292: Mat A,B;
293: PCGetOperators(ksp->B,&A,PETSC_NULL,PETSC_NULL);
294: MatComputeExplicitOperator(A,&B);
295: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(ksp->comm),PETSC_VIEWER_ASCII_MATLAB);
296: MatView(B,PETSC_VIEWER_STDOUT_(ksp->comm));
297: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(ksp->comm));
298: MatDestroy(B);
299: }
300: return(0);
301: }
303: /*@
304: KSPSolveTranspose - Solves the transpose of a linear system. Usually
305: accessed through SLESSolveTranspose().
307: Collective on KSP
309: Input Parameter:
310: . ksp - iterative context obtained from KSPCreate()
312: Output Parameter:
313: . its - number of iterations required
315: Notes:
316: On return, the parameter "its" contains either the iteration
317: number at which convergence was successfully reached, or the
318: negative of the iteration at which divergence or breakdown was detected.
320: Currently only supported by KSPType of KSPPREONLY. This routine is usally
321: only used internally by the BiCG solver on the subblocks in BJacobi and ASM.
323: Level: developer
325: .keywords: KSP, solve, linear system
327: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPDefaultConverged(),
328: SLESSolve(), SLESGetKSP()
329: @*/
330: int KSPSolveTranspose(KSP ksp,int *its)
331: {
332: int ierr;
333: PetscScalar zero = 0.0;
338: if (!ksp->setupcalled){ KSPSetUp(ksp);}
339: if (ksp->guess_zero) { VecSet(&zero,ksp->vec_sol);}
340: ksp->transpose_solve = PETSC_TRUE;
341: (*ksp->ops->solve)(ksp,its);
342: return(0);
343: }
345: /*@C
346: KSPDestroy - Destroys KSP context.
348: Collective on KSP
350: Input Parameter:
351: . ksp - iterative context obtained from KSPCreate()
353: Level: developer
355: .keywords: KSP, destroy
357: .seealso: KSPCreate(), KSPSetUp(), KSPSolve()
358: @*/
359: int KSPDestroy(KSP ksp)
360: {
361: int i,ierr;
365: if (--ksp->refct > 0) return(0);
367: /* if memory was published with AMS then destroy it */
368: PetscObjectDepublish(ksp);
370: if (ksp->ops->destroy) {
371: (*ksp->ops->destroy)(ksp);
372: }
373: for (i=0; i<ksp->numbermonitors; i++) {
374: if (ksp->monitordestroy[i]) {
375: (*ksp->monitordestroy[i])(ksp->monitorcontext[i]);
376: }
377: }
378: PetscLogObjectDestroy(ksp);
379: PetscHeaderDestroy(ksp);
380: return(0);
381: }
383: /*@
384: KSPSetPreconditionerSide - Sets the preconditioning side.
386: Collective on KSP
388: Input Parameter:
389: . ksp - iterative context obtained from KSPCreate()
391: Output Parameter:
392: . side - the preconditioning side, where side is one of
393: .vb
394: PC_LEFT - left preconditioning (default)
395: PC_RIGHT - right preconditioning
396: PC_SYMMETRIC - symmetric preconditioning
397: .ve
399: Options Database Keys:
400: + -ksp_left_pc - Sets left preconditioning
401: . -ksp_right_pc - Sets right preconditioning
402: - -ksp_symmetric_pc - Sets symmetric preconditioning
404: Notes:
405: Left preconditioning is used by default. Symmetric preconditioning is
406: currently available only for the KSPQCG method. Note, however, that
407: symmetric preconditioning can be emulated by using either right or left
408: preconditioning and a pre or post processing step.
410: Level: intermediate
412: .keywords: KSP, set, right, left, symmetric, side, preconditioner, flag
414: .seealso: KSPGetPreconditionerSide()
415: @*/
416: int KSPSetPreconditionerSide(KSP ksp,PCSide side)
417: {
420: ksp->pc_side = side;
421: return(0);
422: }
424: /*@C
425: KSPGetPreconditionerSide - Gets the preconditioning side.
427: Not Collective
429: Input Parameter:
430: . ksp - iterative context obtained from KSPCreate()
432: Output Parameter:
433: . side - the preconditioning side, where side is one of
434: .vb
435: PC_LEFT - left preconditioning (default)
436: PC_RIGHT - right preconditioning
437: PC_SYMMETRIC - symmetric preconditioning
438: .ve
440: Level: intermediate
442: .keywords: KSP, get, right, left, symmetric, side, preconditioner, flag
444: .seealso: KSPSetPreconditionerSide()
445: @*/
446: int KSPGetPreconditionerSide(KSP ksp,PCSide *side)
447: {
450: *side = ksp->pc_side;
451: return(0);
452: }
454: /*@
455: KSPGetTolerances - Gets the relative, absolute, divergence, and maximum
456: iteration tolerances used by the default KSP convergence tests.
458: Not Collective
460: Input Parameter:
461: . ksp - the Krylov subspace context
462:
463: Output Parameters:
464: + rtol - the relative convergence tolerance
465: . atol - the absolute convergence tolerance
466: . dtol - the divergence tolerance
467: - maxits - maximum number of iterations
469: Notes:
470: The user can specify PETSC_NULL for any parameter that is not needed.
472: Level: intermediate
474: .keywords: KSP, get, tolerance, absolute, relative, divergence, convergence,
475: maximum, iterations
477: .seealso: KSPSetTolerances()
478: @*/
479: int KSPGetTolerances(KSP ksp,PetscReal *rtol,PetscReal *atol,PetscReal *dtol,int *maxits)
480: {
483: if (atol) *atol = ksp->atol;
484: if (rtol) *rtol = ksp->rtol;
485: if (dtol) *dtol = ksp->divtol;
486: if (maxits) *maxits = ksp->max_it;
487: return(0);
488: }
490: /*@
491: KSPSetTolerances - Sets the relative, absolute, divergence, and maximum
492: iteration tolerances used by the default KSP convergence testers.
494: Collective on KSP
496: Input Parameters:
497: + ksp - the Krylov subspace context
498: . rtol - the relative convergence tolerance
499: (relative decrease in the residual norm)
500: . atol - the absolute convergence tolerance
501: (absolute size of the residual norm)
502: . dtol - the divergence tolerance
503: (amount residual can increase before KSPDefaultConverged()
504: concludes that the method is diverging)
505: - maxits - maximum number of iterations to use
507: Options Database Keys:
508: + -ksp_atol <atol> - Sets atol
509: . -ksp_rtol <rtol> - Sets rtol
510: . -ksp_divtol <dtol> - Sets dtol
511: - -ksp_max_it <maxits> - Sets maxits
513: Notes:
514: Use PETSC_DEFAULT to retain the default value of any of the tolerances.
516: See KSPDefaultConverged() for details on the use of these parameters
517: in the default convergence test. See also KSPSetConvergenceTest()
518: for setting user-defined stopping criteria.
520: Level: intermediate
522: .keywords: KSP, set, tolerance, absolute, relative, divergence,
523: convergence, maximum, iterations
525: .seealso: KSPGetTolerances(), KSPDefaultConverged(), KSPSetConvergenceTest()
526: @*/
527: int KSPSetTolerances(KSP ksp,PetscReal rtol,PetscReal atol,PetscReal dtol,int maxits)
528: {
531: if (atol != PETSC_DEFAULT) ksp->atol = atol;
532: if (rtol != PETSC_DEFAULT) ksp->rtol = rtol;
533: if (dtol != PETSC_DEFAULT) ksp->divtol = dtol;
534: if (maxits != PETSC_DEFAULT) ksp->max_it = maxits;
535: return(0);
536: }
538: /*@
539: KSPSetInitialGuessNonzero - Tells the iterative solver that the
540: initial guess is nonzero; otherwise KSP assumes the initial guess
541: is to be zero (and thus zeros it out before solving).
543: Collective on KSP
545: Input Parameters:
546: + ksp - iterative context obtained from SLESGetKSP() or KSPCreate()
547: - flg - PETSC_TRUE or PETSC_FALSE
549: Level: beginner
551: Notes:
552: If this is not called the X vector is zeroed in the call to
553: SLESSolve() (or KSPSolve()).
555: .keywords: KSP, set, initial guess, nonzero
557: .seealso: KSPGetIntialGuessNonzero()
558: @*/
559: int KSPSetInitialGuessNonzero(KSP ksp,PetscTruth flg)
560: {
562: ksp->guess_zero = (PetscTruth)!(int)flg;
563: return(0);
564: }
566: /*@
567: KSPGetInitialGuessNonzero - Determines whether the KSP solver is using
568: a zero initial guess.
570: Not Collective
572: Input Parameter:
573: . ksp - iterative context obtained from KSPCreate()
575: Output Parameter:
576: . flag - PETSC_TRUE if guess is nonzero, else PETSC_FALSE
578: Level: intermediate
580: .keywords: KSP, set, initial guess, nonzero
582: .seealso: KSPSetIntialGuessNonzero()
583: @*/
584: int KSPGetInitialGuessNonzero(KSP ksp,PetscTruth *flag)
585: {
587: if (ksp->guess_zero) *flag = PETSC_FALSE;
588: else *flag = PETSC_TRUE;
589: return(0);
590: }
592: /*@
593: KSPSetComputeSingularValues - Sets a flag so that the extreme singular
594: values will be calculated via a Lanczos or Arnoldi process as the linear
595: system is solved.
597: Collective on KSP
599: Input Parameters:
600: + ksp - iterative context obtained from KSPCreate()
601: - flg - PETSC_TRUE or PETSC_FALSE
603: Options Database Key:
604: . -ksp_singmonitor - Activates KSPSetComputeSingularValues()
606: Notes:
607: Currently this option is not valid for all iterative methods.
609: Many users may just want to use the monitoring routine
610: KSPSingularValueMonitor() (which can be set with option -ksp_singmonitor)
611: to print the singular values at each iteration of the linear solve.
613: Level: advanced
615: .keywords: KSP, set, compute, singular values
617: .seealso: KSPComputeExtremeSingularValues(), KSPSingularValueMonitor()
618: @*/
619: int KSPSetComputeSingularValues(KSP ksp,PetscTruth flg)
620: {
623: ksp->calc_sings = flg;
624: return(0);
625: }
627: /*@
628: KSPSetComputeEigenvalues - Sets a flag so that the extreme eigenvalues
629: values will be calculated via a Lanczos or Arnoldi process as the linear
630: system is solved.
632: Collective on KSP
634: Input Parameters:
635: + ksp - iterative context obtained from KSPCreate()
636: - flg - PETSC_TRUE or PETSC_FALSE
638: Notes:
639: Currently this option is not valid for all iterative methods.
641: Level: advanced
643: .keywords: KSP, set, compute, eigenvalues
645: .seealso: KSPComputeEigenvalues(), KSPComputeEigenvaluesExplicitly()
646: @*/
647: int KSPSetComputeEigenvalues(KSP ksp,PetscTruth flg)
648: {
651: ksp->calc_sings = flg;
652: return(0);
653: }
655: /*@
656: KSPSetRhs - Sets the right-hand-side vector for the linear system to
657: be solved.
659: Collective on KSP and Vec
661: Input Parameters:
662: + ksp - iterative context obtained from KSPCreate()
663: - b - right-hand-side vector
665: Level: developer
667: .keywords: KSP, set, right-hand-side, rhs
669: .seealso: KSPGetRhs(), KSPSetSolution()
670: @*/
671: int KSPSetRhs(KSP ksp,Vec b)
672: {
676: ksp->vec_rhs = (b);
677: return(0);
678: }
680: /*@C
681: KSPGetRhs - Gets the right-hand-side vector for the linear system to
682: be solved.
684: Not Collective
686: Input Parameter:
687: . ksp - iterative context obtained from KSPCreate()
689: Output Parameter:
690: . r - right-hand-side vector
692: Level: developer
694: .keywords: KSP, get, right-hand-side, rhs
696: .seealso: KSPSetRhs(), KSPGetSolution()
697: @*/
698: int KSPGetRhs(KSP ksp,Vec *r)
699: {
702: *r = ksp->vec_rhs;
703: return(0);
704: }
706: /*@
707: KSPSetSolution - Sets the location of the solution for the
708: linear system to be solved.
710: Collective on KSP and Vec
712: Input Parameters:
713: + ksp - iterative context obtained from KSPCreate()
714: - x - solution vector
716: Level: developer
718: .keywords: KSP, set, solution
720: .seealso: KSPSetRhs(), KSPGetSolution()
721: @*/
722: int KSPSetSolution(KSP ksp,Vec x)
723: {
727: ksp->vec_sol = (x);
728: return(0);
729: }
731: /*@C
732: KSPGetSolution - Gets the location of the solution for the
733: linear system to be solved. Note that this may not be where the solution
734: is stored during the iterative process; see KSPBuildSolution().
736: Not Collective
738: Input Parameters:
739: . ksp - iterative context obtained from KSPCreate()
741: Output Parameters:
742: . v - solution vector
744: Level: developer
746: .keywords: KSP, get, solution
748: .seealso: KSPGetRhs(), KSPSetSolution(), KSPBuildSolution()
749: @*/
750: int KSPGetSolution(KSP ksp,Vec *v)
751: {
754: *v = ksp->vec_sol;
755: return(0);
756: }
758: /*@
759: KSPSetPC - Sets the preconditioner to be used to calculate the
760: application of the preconditioner on a vector.
762: Collective on KSP
764: Input Parameters:
765: + ksp - iterative context obtained from KSPCreate()
766: - B - the preconditioner object
768: Notes:
769: Use KSPGetPC() to retrieve the preconditioner context (for example,
770: to free it at the end of the computations).
772: Level: developer
774: .keywords: KSP, set, precondition, Binv
776: .seealso: KSPGetPC()
777: @*/
778: int KSPSetPC(KSP ksp,PC B)
779: {
784: ksp->B = B;
785: return(0);
786: }
788: /*@C
789: KSPGetPC - Returns a pointer to the preconditioner context
790: set with KSPSetPC().
792: Not Collective
794: Input Parameters:
795: . ksp - iterative context obtained from KSPCreate()
797: Output Parameter:
798: . B - preconditioner context
800: Level: developer
802: .keywords: KSP, get, preconditioner, Binv
804: .seealso: KSPSetPC()
805: @*/
806: int KSPGetPC(KSP ksp,PC *B)
807: {
810: *B = ksp->B;
811: return(0);
812: }
814: /*@C
815: KSPSetMonitor - Sets an ADDITIONAL function to be called at every iteration to monitor
816: the residual/error etc.
817:
818: Collective on KSP
820: Input Parameters:
821: + ksp - iterative context obtained from KSPCreate()
822: . monitor - pointer to function (if this is PETSC_NULL, it turns off monitoring
823: . mctx - [optional] context for private data for the
824: monitor routine (use PETSC_NULL if no context is desired)
825: - monitordestroy - [optional] routine that frees monitor context
826: (may be PETSC_NULL)
828: Calling Sequence of monitor:
829: $ monitor (KSP ksp, int it, PetscReal rnorm, void *mctx)
831: + ksp - iterative context obtained from KSPCreate()
832: . it - iteration number
833: . rnorm - (estimated) 2-norm of (preconditioned) residual
834: - mctx - optional monitoring context, as set by KSPSetMonitor()
836: Options Database Keys:
837: + -ksp_monitor - sets KSPDefaultMonitor()
838: . -ksp_truemonitor - sets KSPTrueMonitor()
839: . -ksp_xmonitor - sets line graph monitor,
840: uses KSPLGMonitorCreate()
841: . -ksp_xtruemonitor - sets line graph monitor,
842: uses KSPLGMonitorCreate()
843: . -ksp_singmonitor - sets KSPSingularValueMonitor()
844: - -ksp_cancelmonitors - cancels all monitors that have
845: been hardwired into a code by
846: calls to KSPSetMonitor(), but
847: does not cancel those set via
848: the options database.
850: Notes:
851: The default is to do nothing. To print the residual, or preconditioned
852: residual if KSPSetNormType(ksp,KSP_PRECONDITIONED_NORM) was called, use
853: KSPDefaultMonitor() as the monitoring routine, with a null monitoring
854: context.
856: Several different monitoring routines may be set by calling
857: KSPSetMonitor() multiple times; all will be called in the
858: order in which they were set.
860: Level: beginner
862: .keywords: KSP, set, monitor
864: .seealso: KSPDefaultMonitor(), KSPLGMonitorCreate(), KSPClearMonitor()
865: @*/
866: int KSPSetMonitor(KSP ksp,int (*monitor)(KSP,int,PetscReal,void*),void *mctx,int (*monitordestroy)(void*))
867: {
870: if (ksp->numbermonitors >= MAXKSPMONITORS) {
871: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many KSP monitors set");
872: }
873: ksp->monitor[ksp->numbermonitors] = monitor;
874: ksp->monitordestroy[ksp->numbermonitors] = monitordestroy;
875: ksp->monitorcontext[ksp->numbermonitors++] = (void*)mctx;
876: return(0);
877: }
879: /*@
880: KSPClearMonitor - Clears all monitors for a KSP object.
882: Collective on KSP
884: Input Parameters:
885: . ksp - iterative context obtained from KSPCreate()
887: Options Database Key:
888: . -ksp_cancelmonitors - Cancels all monitors that have
889: been hardwired into a code by calls to KSPSetMonitor(),
890: but does not cancel those set via the options database.
892: Level: intermediate
894: .keywords: KSP, set, monitor
896: .seealso: KSPDefaultMonitor(), KSPLGMonitorCreate(), KSPSetMonitor()
897: @*/
898: int KSPClearMonitor(KSP ksp)
899: {
902: ksp->numbermonitors = 0;
903: return(0);
904: }
906: /*@C
907: KSPGetMonitorContext - Gets the monitoring context, as set by
908: KSPSetMonitor() for the FIRST monitor only.
910: Not Collective
912: Input Parameter:
913: . ksp - iterative context obtained from KSPCreate()
915: Output Parameter:
916: . ctx - monitoring context
918: Level: intermediate
920: .keywords: KSP, get, monitor, context
922: .seealso: KSPDefaultMonitor(), KSPLGMonitorCreate()
923: @*/
924: int KSPGetMonitorContext(KSP ksp,void **ctx)
925: {
928: *ctx = (ksp->monitorcontext[0]);
929: return(0);
930: }
932: /*@
933: KSPSetResidualHistory - Sets the array used to hold the residual history.
934: If set, this array will contain the residual norms computed at each
935: iteration of the solver.
937: Not Collective
939: Input Parameters:
940: + ksp - iterative context obtained from KSPCreate()
941: . a - array to hold history
942: . na - size of a
943: - reset - PETSC_TRUE indicates the history counter is reset to zero
944: for each new linear solve
946: Level: advanced
948: .keywords: KSP, set, residual, history, norm
950: .seealso: KSPGetResidualHistory()
952: @*/
953: int KSPSetResidualHistory(KSP ksp,PetscReal *a,int na,PetscTruth reset)
954: {
959: if (na != PETSC_DECIDE && a != PETSC_NULL) {
960: ksp->res_hist = a;
961: ksp->res_hist_max = na;
962: } else {
963: ksp->res_hist_max = 1000;
964: PetscMalloc(ksp->res_hist_max*sizeof(PetscReal),&ksp->res_hist);
965: }
966: ksp->res_hist_len = 0;
967: ksp->res_hist_reset = reset;
970: return(0);
971: }
973: /*@C
974: KSPGetResidualHistory - Gets the array used to hold the residual history
975: and the number of residuals it contains.
977: Not Collective
979: Input Parameter:
980: . ksp - iterative context obtained from KSPCreate()
982: Output Parameters:
983: + a - pointer to array to hold history (or PETSC_NULL)
984: - na - number of used entries in a (or PETSC_NULL)
986: Level: advanced
988: Notes:
989: Can only call after a KSPSetResidualHistory() otherwise returns 0.
991: The Fortran version of this routine has a calling sequence
992: $ call KSPGetResidualHistory(KSP ksp, integer na, integer ierr)
994: .keywords: KSP, get, residual, history, norm
996: .seealso: KSPGetResidualHistory()
998: @*/
999: int KSPGetResidualHistory(KSP ksp,PetscReal **a,int *na)
1000: {
1003: if (a) *a = ksp->res_hist;
1004: if (na) *na = ksp->res_hist_len;
1005: return(0);
1006: }
1008: /*@C
1009: KSPSetConvergenceTest - Sets the function to be used to determine
1010: convergence.
1012: Collective on KSP
1014: Input Parameters:
1015: + ksp - iterative context obtained from KSPCreate()
1016: . converge - pointer to int function
1017: - cctx - context for private data for the convergence routine (may be null)
1019: Calling sequence of converge:
1020: $ converge (KSP ksp, int it, PetscReal rnorm, KSPConvergedReason *reason,void *mctx)
1022: + ksp - iterative context obtained from KSPCreate()
1023: . it - iteration number
1024: . rnorm - (estimated) 2-norm of (preconditioned) residual
1025: . reason - the reason why it has converged or diverged
1026: - cctx - optional convergence context, as set by KSPSetConvergenceTest()
1029: Notes:
1030: Must be called after the KSP type has been set so put this after
1031: a call to KSPSetType(), or KSPSetFromOptions().
1033: The default convergence test, KSPDefaultConverged(), aborts if the
1034: residual grows to more than 10000 times the initial residual.
1036: The default is a combination of relative and absolute tolerances.
1037: The residual value that is tested may be an approximation; routines
1038: that need exact values should compute them.
1040: Level: advanced
1042: .keywords: KSP, set, convergence, test, context
1044: .seealso: KSPDefaultConverged(), KSPGetConvergenceContext()
1045: @*/
1046: int KSPSetConvergenceTest(KSP ksp,int (*converge)(KSP,int,PetscReal,KSPConvergedReason*,void*),void *cctx)
1047: {
1050: ksp->converged = converge;
1051: ksp->cnvP = (void*)cctx;
1052: return(0);
1053: }
1055: /*@C
1056: KSPGetConvergenceContext - Gets the convergence context set with
1057: KSPSetConvergenceTest().
1059: Not Collective
1061: Input Parameter:
1062: . ksp - iterative context obtained from KSPCreate()
1064: Output Parameter:
1065: . ctx - monitoring context
1067: Level: advanced
1069: .keywords: KSP, get, convergence, test, context
1071: .seealso: KSPDefaultConverged(), KSPSetConvergenceTest()
1072: @*/
1073: int KSPGetConvergenceContext(KSP ksp,void **ctx)
1074: {
1077: *ctx = ksp->cnvP;
1078: return(0);
1079: }
1081: /*@C
1082: KSPBuildSolution - Builds the approximate solution in a vector provided.
1083: This routine is NOT commonly needed (see SLESSolve()).
1085: Collective on KSP
1087: Input Parameter:
1088: . ctx - iterative context obtained from KSPCreate()
1090: Output Parameter:
1091: Provide exactly one of
1092: + v - location to stash solution.
1093: - V - the solution is returned in this location. This vector is created
1094: internally. This vector should NOT be destroyed by the user with
1095: VecDestroy().
1097: Notes:
1098: This routine can be used in one of two ways
1099: .vb
1100: KSPBuildSolution(ksp,PETSC_NULL,&V);
1101: or
1102: KSPBuildSolution(ksp,v,PETSC_NULL);
1103: .ve
1104: In the first case an internal vector is allocated to store the solution
1105: (the user cannot destroy this vector). In the second case the solution
1106: is generated in the vector that the user provides. Note that for certain
1107: methods, such as KSPCG, the second case requires a copy of the solution,
1108: while in the first case the call is essentially free since it simply
1109: returns the vector where the solution already is stored.
1111: Level: advanced
1113: .keywords: KSP, build, solution
1115: .seealso: KSPGetSolution(), KSPBuildResidual()
1116: @*/
1117: int KSPBuildSolution(KSP ksp,Vec v,Vec *V)
1118: {
1123: if (!V && !v) SETERRQ(PETSC_ERR_ARG_WRONG,"Must provide either v or V");
1124: if (!V) V = &v;
1125: (*ksp->ops->buildsolution)(ksp,v,V);
1126: return(0);
1127: }
1129: /*@C
1130: KSPBuildResidual - Builds the residual in a vector provided.
1132: Collective on KSP
1134: Input Parameter:
1135: . ksp - iterative context obtained from KSPCreate()
1137: Output Parameters:
1138: + v - optional location to stash residual. If v is not provided,
1139: then a location is generated.
1140: . t - work vector. If not provided then one is generated.
1141: - V - the residual
1143: Notes:
1144: Regardless of whether or not v is provided, the residual is
1145: returned in V.
1147: Level: advanced
1149: .keywords: KSP, build, residual
1151: .seealso: KSPBuildSolution()
1152: @*/
1153: int KSPBuildResidual(KSP ksp,Vec t,Vec v,Vec *V)
1154: {
1155: int flag = 0,ierr;
1156: Vec w = v,tt = t;
1160: if (!w) {
1161: VecDuplicate(ksp->vec_rhs,&w);
1162: PetscLogObjectParent((PetscObject)ksp,w);
1163: }
1164: if (!tt) {
1165: VecDuplicate(ksp->vec_rhs,&tt); flag = 1;
1166: PetscLogObjectParent((PetscObject)ksp,tt);
1167: }
1168: (*ksp->ops->buildresidual)(ksp,tt,w,V);
1169: if (flag) {VecDestroy(tt);}
1170: return(0);
1171: }