Actual source code: iterativ.c
petsc-3.4.5 2014-06-29
1: /*
2: This file contains some simple default routines.
3: These routines should be SHORT, since they will be included in every
4: executable image that uses the iterative routines (note that, through
5: the registry system, we provide a way to load only the truely necessary
6: files)
7: */
8: #include <petsc-private/kspimpl.h> /*I "petscksp.h" I*/
9: #include <petscdmshell.h>
13: /*@
14: KSPGetResidualNorm - Gets the last (approximate preconditioned)
15: residual norm that has been computed.
17: Not Collective
19: Input Parameters:
20: . ksp - the iterative context
22: Output Parameters:
23: . rnorm - residual norm
25: Level: intermediate
27: .keywords: KSP, get, residual norm
29: .seealso: KSPBuildResidual()
30: @*/
31: PetscErrorCode KSPGetResidualNorm(KSP ksp,PetscReal *rnorm)
32: {
36: *rnorm = ksp->rnorm;
37: return(0);
38: }
42: /*@
43: KSPGetIterationNumber - Gets the current iteration number; if the
44: KSPSolve() is complete, returns the number of iterations
45: used.
47: Not Collective
49: Input Parameters:
50: . ksp - the iterative context
52: Output Parameters:
53: . its - number of iterations
55: Level: intermediate
57: Notes:
58: During the ith iteration this returns i-1
59: .keywords: KSP, get, residual norm
61: .seealso: KSPBuildResidual(), KSPGetResidualNorm()
62: @*/
63: PetscErrorCode KSPGetIterationNumber(KSP ksp,PetscInt *its)
64: {
68: *its = ksp->its;
69: return(0);
70: }
74: /*@C
75: KSPMonitorSingularValue - Prints the two norm of the true residual and
76: estimation of the extreme singular values of the preconditioned problem
77: at each iteration.
79: Logically Collective on KSP
81: Input Parameters:
82: + ksp - the iterative context
83: . n - the iteration
84: - rnorm - the two norm of the residual
86: Options Database Key:
87: . -ksp_monitor_singular_value - Activates KSPMonitorSingularValue()
89: Notes:
90: The CG solver uses the Lanczos technique for eigenvalue computation,
91: while GMRES uses the Arnoldi technique; other iterative methods do
92: not currently compute singular values.
94: Level: intermediate
96: .keywords: KSP, CG, default, monitor, extreme, singular values, Lanczos, Arnoldi
98: .seealso: KSPComputeExtremeSingularValues()
99: @*/
100: PetscErrorCode KSPMonitorSingularValue(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy)
101: {
102: PetscReal emin,emax,c;
104: PetscViewer viewer = (PetscViewer) dummy;
108: if (!viewer) {
109: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ksp),&viewer);
110: }
111: PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
112: if (!ksp->calc_sings) {
113: PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm %14.12e \n",n,(double)rnorm);
114: } else {
115: KSPComputeExtremeSingularValues(ksp,&emax,&emin);
116: c = emax/emin;
117: PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm %14.12e %% max %14.12e min %14.12e max/min %14.12e\n",n,(double)rnorm,(double)emax,(double)emin,(double)c);
118: }
119: PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
120: return(0);
121: }
125: /*@C
126: KSPMonitorSolution - Monitors progress of the KSP solvers by calling
127: VecView() for the approximate solution at each iteration.
129: Collective on KSP
131: Input Parameters:
132: + ksp - the KSP context
133: . its - iteration number
134: . fgnorm - 2-norm of residual (or gradient)
135: - dummy - either a viewer or NULL
137: Level: intermediate
139: Notes:
140: For some Krylov methods such as GMRES constructing the solution at
141: each iteration is expensive, hence using this will slow the code.
143: .keywords: KSP, nonlinear, vector, monitor, view
145: .seealso: KSPMonitorSet(), KSPMonitorDefault(), VecView()
146: @*/
147: PetscErrorCode KSPMonitorSolution(KSP ksp,PetscInt its,PetscReal fgnorm,void *dummy)
148: {
150: Vec x;
151: PetscViewer viewer = (PetscViewer) dummy;
154: KSPBuildSolution(ksp,NULL,&x);
155: if (!viewer) {
156: MPI_Comm comm;
157: PetscObjectGetComm((PetscObject)ksp,&comm);
158: viewer = PETSC_VIEWER_DRAW_(comm);
159: }
160: VecView(x,viewer);
161: return(0);
162: }
166: /*@C
167: KSPMonitorDefault - Print the residual norm at each iteration of an
168: iterative solver.
170: Collective on KSP
172: Input Parameters:
173: + ksp - iterative context
174: . n - iteration number
175: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
176: - dummy - unused monitor context
178: Level: intermediate
180: .keywords: KSP, default, monitor, residual
182: .seealso: KSPMonitorSet(), KSPMonitorTrueResidualNorm(), KSPMonitorLGResidualNormCreate()
183: @*/
184: PetscErrorCode KSPMonitorDefault(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy)
185: {
187: PetscViewer viewer = (PetscViewer) dummy;
190: if (!viewer) {
191: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ksp),&viewer);
192: }
193: PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
194: if (n == 0 && ((PetscObject)ksp)->prefix) {
195: PetscViewerASCIIPrintf(viewer," Residual norms for %s solve.\n",((PetscObject)ksp)->prefix);
196: }
197: PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm %14.12e \n",n,(double)rnorm);
198: PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
199: return(0);
200: }
204: /*@C
205: KSPMonitorTrueResidualNorm - Prints the true residual norm as well as the preconditioned
206: residual norm at each iteration of an iterative solver.
208: Collective on KSP
210: Input Parameters:
211: + ksp - iterative context
212: . n - iteration number
213: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
214: - dummy - unused monitor context
216: Options Database Key:
217: . -ksp_monitor_true_residual - Activates KSPMonitorTrueResidualNorm()
219: Notes:
220: When using right preconditioning, these values are equivalent.
222: Level: intermediate
224: .keywords: KSP, default, monitor, residual
226: .seealso: KSPMonitorSet(), KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(),KSPMonitorTrueResidualMaxNorm()
227: @*/
228: PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy)
229: {
231: Vec resid;
232: PetscReal truenorm,bnorm;
233: PetscViewer viewer = (PetscViewer)dummy;
234: char normtype[256];
237: if (!viewer) {
238: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ksp),&viewer);
239: }
240: PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
241: if (n == 0 && ((PetscObject)ksp)->prefix) {
242: PetscViewerASCIIPrintf(viewer," Residual norms for %s solve.\n",((PetscObject)ksp)->prefix);
243: }
244: KSPBuildResidual(ksp,NULL,NULL,&resid);
245: VecNorm(resid,NORM_2,&truenorm);
246: VecDestroy(&resid);
247: VecNorm(ksp->vec_rhs,NORM_2,&bnorm);
248: PetscStrncpy(normtype,KSPNormTypes[ksp->normtype],sizeof(normtype));
249: PetscStrtolower(normtype);
250: PetscViewerASCIIPrintf(viewer,"%3D KSP %s resid norm %14.12e true resid norm %14.12e ||r(i)||/||b|| %14.12e\n",n,normtype,(double)rnorm,(double)truenorm,(double)(truenorm/bnorm));
251: PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
252: return(0);
253: }
257: /*@C
258: KSPMonitorTrueResidualMaxNorm - Prints the true residual max norm as well as the preconditioned
259: residual norm at each iteration of an iterative solver.
261: Collective on KSP
263: Input Parameters:
264: + ksp - iterative context
265: . n - iteration number
266: . rnorm - norm (preconditioned) residual value (may be estimated).
267: - dummy - unused monitor context
269: Options Database Key:
270: . -ksp_monitor_max - Activates KSPMonitorTrueResidualMaxNorm()
272: Notes:
273: This could be implemented (better) with a flag in ksp.
275: Level: intermediate
277: .keywords: KSP, default, monitor, residual
279: .seealso: KSPMonitorSet(), KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(),KSPMonitorTrueResidualNorm()
280: @*/
281: PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy)
282: {
284: Vec resid;
285: PetscReal truenorm,bnorm;
286: PetscViewer viewer = (PetscViewer)dummy;
287: char normtype[256];
290: if (!viewer) {
291: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ksp),&viewer);
292: }
293: PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
294: if (n == 0 && ((PetscObject)ksp)->prefix) {
295: PetscViewerASCIIPrintf(viewer," Residual norms (max) for %s solve.\n",((PetscObject)ksp)->prefix);
296: }
297: KSPBuildResidual(ksp,NULL,NULL,&resid);
298: VecNorm(resid,NORM_INFINITY,&truenorm);
299: VecDestroy(&resid);
300: VecNorm(ksp->vec_rhs,NORM_INFINITY,&bnorm);
301: PetscStrncpy(normtype,KSPNormTypes[ksp->normtype],sizeof(normtype));
302: PetscStrtolower(normtype);
303: /* PetscViewerASCIIPrintf(viewer,"%3D KSP %s resid norm %14.12e true resid norm %14.12e ||r(i)||_inf/||b||_inf %14.12e\n",n,normtype,(double)rnorm,(double)truenorm,(double)(truenorm/bnorm)); */
304: PetscViewerASCIIPrintf(viewer,"%3D KSP true resid max norm %14.12e ||r(i)||/||b|| %14.12e\n",n,(double)truenorm,(double)(truenorm/bnorm));
305: PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
306: return(0);
307: }
311: PetscErrorCode KSPMonitorRange_Private(KSP ksp,PetscInt it,PetscReal *per)
312: {
314: Vec resid;
315: PetscReal rmax,pwork;
316: PetscInt i,n,N;
317: const PetscScalar *r;
320: KSPBuildResidual(ksp,NULL,NULL,&resid);
321: VecNorm(resid,NORM_INFINITY,&rmax);
322: VecGetLocalSize(resid,&n);
323: VecGetSize(resid,&N);
324: VecGetArrayRead(resid,&r);
325: pwork = 0.0;
326: for (i=0; i<n; i++) pwork += (PetscAbsScalar(r[i]) > .20*rmax);
327: VecRestoreArrayRead(resid,&r);
328: VecDestroy(&resid);
329: MPI_Allreduce(&pwork,per,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ksp));
330: *per = *per/N;
331: return(0);
332: }
336: /*@C
337: KSPMonitorRange - Prints the percentage of residual elements that are more then 10 percent of the maximum value.
339: Collective on KSP
341: Input Parameters:
342: + ksp - iterative context
343: . it - iteration number
344: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
345: - dummy - unused monitor context
347: Options Database Key:
348: . -ksp_monitor_range - Activates KSPMonitorRange()
350: Level: intermediate
352: .keywords: KSP, default, monitor, residual
354: .seealso: KSPMonitorSet(), KSPMonitorDefault(), KSPMonitorLGResidualNormCreate()
355: @*/
356: PetscErrorCode KSPMonitorRange(KSP ksp,PetscInt it,PetscReal rnorm,void *dummy)
357: {
358: PetscErrorCode ierr;
359: PetscReal perc,rel;
360: PetscViewer viewer = (PetscViewer)dummy;
361: /* should be in a MonitorRangeContext */
362: static PetscReal prev;
365: if (!viewer) {
366: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ksp),&viewer);
367: }
368: PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
369: if (!it) prev = rnorm;
370: if (it == 0 && ((PetscObject)ksp)->prefix) {
371: PetscViewerASCIIPrintf(viewer," Residual norms for %s solve.\n",((PetscObject)ksp)->prefix);
372: }
373: KSPMonitorRange_Private(ksp,it,&perc);
375: rel = (prev - rnorm)/prev;
376: prev = rnorm;
377: PetscViewerASCIIPrintf(viewer,"%3D KSP preconditioned resid norm %14.12e Percent values above 20 percent of maximum %5.2f relative decrease %5.2e ratio %5.2e \n",it,(double)rnorm,(double)100.0*perc,(double)rel,(double)rel/perc);
378: PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
379: return(0);
380: }
384: /*@C
385: KSPMonitorDynamicTolerance - Recompute the inner tolerance in every
386: outer iteration in an adaptive way.
388: Collective on KSP
390: Input Parameters:
391: + ksp - iterative context
392: . n - iteration number (not used)
393: . fnorm - the current residual norm
394: . dummy - some context as a C struct. fields:
395: coef: a scaling coefficient. default 1.0. can be passed through
396: -sub_ksp_dynamic_tolerance_param
397: bnrm: norm of the right-hand side. store it to avoid repeated calculation
399: Notes:
400: This may be useful for a flexibly preconditioner Krylov method to
401: control the accuracy of the inner solves needed to gaurantee the
402: convergence of the outer iterations.
404: Level: advanced
406: .keywords: KSP, inner tolerance
408: .seealso: KSPMonitorDynamicToleranceDestroy()
409: @*/
410: PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy)
411: {
413: PC pc;
414: PetscReal outer_rtol, outer_abstol, outer_dtol, inner_rtol;
415: PetscInt outer_maxits,nksp,first,i;
416: KSPDynTolCtx *scale = (KSPDynTolCtx*)dummy;
417: KSP kspinner = NULL, *subksp = NULL;
420: KSPGetPC(ksp, &pc);
422: /* compute inner_rtol */
423: if (scale->bnrm < 0.0) {
424: Vec b;
425: KSPGetRhs(ksp, &b);
426: VecNorm(b, NORM_2, &(scale->bnrm));
427: }
428: KSPGetTolerances(ksp, &outer_rtol, &outer_abstol, &outer_dtol, &outer_maxits);
429: inner_rtol = PetscMin(scale->coef * scale->bnrm * outer_rtol / fnorm, 0.999);
430: /*PetscPrintf(PETSC_COMM_WORLD, " Inner rtol = %G\n", inner_rtol);*/
432: /* if pc is ksp */
433: PCKSPGetKSP(pc, &kspinner);
434: if (kspinner) {
435: KSPSetTolerances(kspinner, inner_rtol, outer_abstol, outer_dtol, outer_maxits);
436: return(0);
437: }
439: /* if pc is bjacobi */
440: PCBJacobiGetSubKSP(pc, &nksp, &first, &subksp);
441: if (subksp) {
442: for (i=0; i<nksp; i++) {
443: KSPSetTolerances(subksp[i], inner_rtol, outer_abstol, outer_dtol, outer_maxits);
444: }
445: return(0);
446: }
448: /* todo: dynamic tolerance may apply to other types of pc too */
449: return(0);
450: }
454: /*
455: Destroy the dummy context used in KSPMonitorDynamicTolerance()
456: */
457: PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy)
458: {
462: PetscFree(*dummy);
463: return(0);
464: }
468: /*
469: Default (short) KSP Monitor, same as KSPMonitorDefault() except
470: it prints fewer digits of the residual as the residual gets smaller.
471: This is because the later digits are meaningless and are often
472: different on different machines; by using this routine different
473: machines will usually generate the same output.
474: */
475: PetscErrorCode KSPMonitorDefaultShort(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy)
476: {
478: PetscViewer viewer = (PetscViewer)dummy;
481: if (!viewer) {
482: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ksp),&viewer);
483: }
484: PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
485: if (its == 0 && ((PetscObject)ksp)->prefix) {
486: PetscViewerASCIIPrintf(viewer," Residual norms for %s solve.\n",((PetscObject)ksp)->prefix);
487: }
489: if (fnorm > 1.e-9) {
490: PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm %G \n",its,fnorm);
491: } else if (fnorm > 1.e-11) {
492: PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm %5.3e \n",its,(double)fnorm);
493: } else {
494: PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm < 1.e-11\n",its);
495: }
496: PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
497: return(0);
498: }
502: /*@C
503: KSPSkipConverged - Convergence test that do not return as converged
504: until the maximum number of iterations is reached.
506: Collective on KSP
508: Input Parameters:
509: + ksp - iterative context
510: . n - iteration number
511: . rnorm - 2-norm residual value (may be estimated)
512: - dummy - unused convergence context
514: Returns:
515: . reason - KSP_CONVERGED_ITERATING, KSP_CONVERGED_ITS
517: Notes:
518: This should be used as the convergence test with the option
519: KSPSetNormType(ksp,KSP_NORM_NONE), since norms of the residual are
520: not computed. Convergence is then declared after the maximum number
521: of iterations have been reached. Useful when one is using CG or
522: BiCGStab as a smoother.
524: Level: advanced
526: .keywords: KSP, default, convergence, residual
528: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPSetNormType()
529: @*/
530: PetscErrorCode KSPSkipConverged(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *dummy)
531: {
535: *reason = KSP_CONVERGED_ITERATING;
536: if (n >= ksp->max_it) *reason = KSP_CONVERGED_ITS;
537: return(0);
538: }
543: /*@C
544: KSPDefaultConvergedCreate - Creates and initializes the space used by the KSPDefaultConverged() function context
546: Collective on KSP
548: Output Parameter:
549: . ctx - convergence context
551: Level: intermediate
553: .keywords: KSP, default, convergence, residual
555: .seealso: KSPDefaultConverged(), KSPDefaultConvergedDestroy(), KSPSetConvergenceTest(), KSPSetTolerances(),
556: KSPSkipConverged(), KSPConvergedReason, KSPGetConvergedReason(), KSPDefaultConvergedSetUIRNorm(), KSPDefaultConvergedSetUMIRNorm()
557: @*/
558: PetscErrorCode KSPDefaultConvergedCreate(void **ctx)
559: {
560: PetscErrorCode ierr;
561: KSPDefaultConvergedCtx *cctx;
564: PetscNew(KSPDefaultConvergedCtx,&cctx);
565: *ctx = cctx;
566: return(0);
567: }
571: /*@
572: KSPDefaultConvergedSetUIRNorm - makes the default convergence test use || B*(b - A*(initial guess))||
573: instead of || B*b ||. In the case of right preconditioner or if KSPSetNormType(ksp,KSP_NORM_UNPRECONDIITONED)
574: is used there is no B in the above formula. UIRNorm is short for Use Initial Residual Norm.
576: Collective on KSP
578: Input Parameters:
579: . ksp - iterative context
581: Options Database:
582: . -ksp_converged_use_initial_residual_norm
584: Notes:
585: Use KSPSetTolerances() to alter the defaults for rtol, abstol, dtol.
587: The precise values of reason are macros such as KSP_CONVERGED_RTOL, which
588: are defined in petscksp.h.
590: If the convergence test is not KSPDefaultConverged() then this is ignored.
592: If right preconditioning is being used then B does not appear in the above formula.
595: Level: intermediate
597: .keywords: KSP, default, convergence, residual
599: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPSkipConverged(), KSPConvergedReason, KSPGetConvergedReason(), KSPDefaultConvergedSetUMIRNorm()
600: @*/
601: PetscErrorCode KSPDefaultConvergedSetUIRNorm(KSP ksp)
602: {
603: KSPDefaultConvergedCtx *ctx = (KSPDefaultConvergedCtx*) ksp->cnvP;
607: if (ksp->converged != KSPDefaultConverged) return(0);
608: if (ctx->mininitialrtol) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONGSTATE,"Cannot use KSPDefaultConvergedSetUIRNorm() and KSPDefaultConvergedSetUMIRNorm() together");
609: ctx->initialrtol = PETSC_TRUE;
610: return(0);
611: }
615: /*@
616: KSPDefaultConvergedSetUMIRNorm - makes the default convergence test use min(|| B*(b - A*(initial guess))||,|| B*b ||)
617: In the case of right preconditioner or if KSPSetNormType(ksp,KSP_NORM_UNPRECONDIITONED)
618: is used there is no B in the above formula. UMIRNorm is short for Use Minimum Initial Residual Norm.
620: Collective on KSP
622: Input Parameters:
623: . ksp - iterative context
625: Options Database:
626: . -ksp_converged_use_min_initial_residual_norm
628: Use KSPSetTolerances() to alter the defaults for rtol, abstol, dtol.
630: The precise values of reason are macros such as KSP_CONVERGED_RTOL, which
631: are defined in petscksp.h.
633: Level: intermediate
635: .keywords: KSP, default, convergence, residual
637: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPSkipConverged(), KSPConvergedReason, KSPGetConvergedReason(), KSPDefaultConvergedSetUIRNorm()
638: @*/
639: PetscErrorCode KSPDefaultConvergedSetUMIRNorm(KSP ksp)
640: {
641: KSPDefaultConvergedCtx *ctx = (KSPDefaultConvergedCtx*) ksp->cnvP;
645: if (ksp->converged != KSPDefaultConverged) return(0);
646: if (ctx->initialrtol) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONGSTATE,"Cannot use KSPDefaultConvergedSetUIRNorm() and KSPDefaultConvergedSetUMIRNorm() together");
647: ctx->mininitialrtol = PETSC_TRUE;
648: return(0);
649: }
653: /*@C
654: KSPDefaultConverged - Determines convergence of
655: the iterative solvers (default code).
657: Collective on KSP
659: Input Parameters:
660: + ksp - iterative context
661: . n - iteration number
662: . rnorm - 2-norm residual value (may be estimated)
663: - ctx - convergence context which must be created by KSPDefaultConvergedCreate()
665: reason is set to:
666: + positive - if the iteration has converged;
667: . negative - if residual norm exceeds divergence threshold;
668: - 0 - otherwise.
670: Notes:
671: KSPDefaultConverged() reaches convergence when
672: $ rnorm < MAX (rtol * rnorm_0, abstol);
673: Divergence is detected if
674: $ rnorm > dtol * rnorm_0,
676: where
677: + rtol = relative tolerance,
678: . abstol = absolute tolerance.
679: . dtol = divergence tolerance,
680: - rnorm_0 is the two norm of the right hand side. When initial guess is non-zero you
681: can call KSPDefaultConvergedSetUIRNorm() to use the norm of (b - A*(initial guess))
682: as the starting point for relative norm convergence testing.
684: Use KSPSetTolerances() to alter the defaults for rtol, abstol, dtol.
686: The precise values of reason are macros such as KSP_CONVERGED_RTOL, which
687: are defined in petscksp.h.
689: Level: intermediate
691: .keywords: KSP, default, convergence, residual
693: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPSkipConverged(), KSPConvergedReason, KSPGetConvergedReason(),
694: KSPDefaultConvergedSetUIRNorm(), KSPDefaultConvergedSetUMIRNorm(), KSPDefaultConvergedCreate(), KSPDefaultConvergedDestroy()
695: @*/
696: PetscErrorCode KSPDefaultConverged(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *ctx)
697: {
698: PetscErrorCode ierr;
699: KSPDefaultConvergedCtx *cctx = (KSPDefaultConvergedCtx*) ctx;
700: KSPNormType normtype;
705: *reason = KSP_CONVERGED_ITERATING;
707: KSPGetNormType(ksp,&normtype);
708: if (normtype == KSP_NORM_NONE) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONGSTATE,"Use KSPSkipConverged() with KSPNormType of KSP_NORM_NONE");
710: if (!cctx) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_NULL,"Convergence context must have been created with KSPDefaultConvergedCreate()");
711: if (!n) {
712: /* if user gives initial guess need to compute norm of b */
713: if (!ksp->guess_zero && !cctx->initialrtol) {
714: PetscReal snorm;
715: if (ksp->normtype == KSP_NORM_UNPRECONDITIONED || ksp->pc_side == PC_RIGHT) {
716: PetscInfo(ksp,"user has provided nonzero initial guess, computing 2-norm of RHS\n");
717: VecNorm(ksp->vec_rhs,NORM_2,&snorm); /* <- b'*b */
718: } else {
719: Vec z;
720: /* Should avoid allocating the z vector each time but cannot stash it in cctx because if KSPReset() is called the vector size might change */
721: VecDuplicate(ksp->vec_rhs,&z);
722: KSP_PCApply(ksp,ksp->vec_rhs,z);
723: if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
724: PetscInfo(ksp,"user has provided nonzero initial guess, computing 2-norm of preconditioned RHS\n");
725: VecNorm(z,NORM_2,&snorm); /* dp <- b'*B'*B*b */
726: } else if (ksp->normtype == KSP_NORM_NATURAL) {
727: PetscScalar norm;
728: PetscInfo(ksp,"user has provided nonzero initial guess, computing natural norm of RHS\n");
729: VecDot(ksp->vec_rhs,z,&norm);
730: snorm = PetscSqrtReal(PetscAbsScalar(norm)); /* dp <- b'*B*b */
731: }
732: VecDestroy(&z);
733: }
734: /* handle special case of zero RHS and nonzero guess */
735: if (!snorm) {
736: PetscInfo(ksp,"Special case, user has provided nonzero initial guess and zero RHS\n");
737: snorm = rnorm;
738: }
739: if (cctx->mininitialrtol) ksp->rnorm0 = PetscMin(snorm,rnorm);
740: else ksp->rnorm0 = snorm;
741: } else {
742: ksp->rnorm0 = rnorm;
743: }
744: ksp->ttol = PetscMax(ksp->rtol*ksp->rnorm0,ksp->abstol);
745: }
747: if (n <= ksp->chknorm) return(0);
749: if (PetscIsInfOrNanScalar(rnorm)) {
750: PetscInfo(ksp,"Linear solver has created a not a number (NaN) as the residual norm, declaring divergence \n");
751: *reason = KSP_DIVERGED_NANORINF;
752: } else if (rnorm <= ksp->ttol) {
753: if (rnorm < ksp->abstol) {
754: PetscInfo3(ksp,"Linear solver has converged. Residual norm %14.12e is less than absolute tolerance %14.12e at iteration %D\n",(double)rnorm,(double)ksp->abstol,n);
755: *reason = KSP_CONVERGED_ATOL;
756: } else {
757: if (cctx->initialrtol) {
758: PetscInfo4(ksp,"Linear solver has converged. Residual norm %14.12e is less than relative tolerance %14.12e times initial residual norm %14.12e at iteration %D\n",(double)rnorm,(double)ksp->rtol,(double)ksp->rnorm0,n);
759: } else {
760: PetscInfo4(ksp,"Linear solver has converged. Residual norm %14.12e is less than relative tolerance %14.12e times initial right hand side norm %14.12e at iteration %D\n",(double)rnorm,(double)ksp->rtol,(double)ksp->rnorm0,n);
761: }
762: *reason = KSP_CONVERGED_RTOL;
763: }
764: } else if (rnorm >= ksp->divtol*ksp->rnorm0) {
765: PetscInfo3(ksp,"Linear solver is diverging. Initial right hand size norm %14.12e, current residual norm %14.12e at iteration %D\n",(double)ksp->rnorm0,(double)rnorm,n);
766: *reason = KSP_DIVERGED_DTOL;
767: }
768: return(0);
769: }
773: /*@C
774: KSPDefaultConvergedDestroy - Frees the space used by the KSPDefaultConverged() function context
776: Collective on KSP
778: Input Parameters:
779: . ctx - convergence context
781: Level: intermediate
783: .keywords: KSP, default, convergence, residual
785: .seealso: KSPDefaultConverged(), KSPDefaultConvergedCreate(), KSPSetConvergenceTest(), KSPSetTolerances(), KSPSkipConverged(),
786: KSPConvergedReason, KSPGetConvergedReason(), KSPDefaultConvergedSetUIRNorm(), KSPDefaultConvergedSetUMIRNorm()
787: @*/
788: PetscErrorCode KSPDefaultConvergedDestroy(void *ctx)
789: {
790: PetscErrorCode ierr;
791: KSPDefaultConvergedCtx *cctx = (KSPDefaultConvergedCtx*) ctx;
794: VecDestroy(&cctx->work);
795: PetscFree(ctx);
796: return(0);
797: }
801: /*
802: KSPBuildSolutionDefault - Default code to create/move the solution.
804: Input Parameters:
805: + ksp - iterative context
806: - v - pointer to the user's vector
808: Output Parameter:
809: . V - pointer to a vector containing the solution
811: Level: advanced
813: Developers Note: This is PETSC_EXTERN because it may be used by user written plugin KSP implementations
815: .keywords: KSP, build, solution, default
817: .seealso: KSPGetSolution(), KSPBuildResidualDefault()
818: */
819: PetscErrorCode KSPBuildSolutionDefault(KSP ksp,Vec v,Vec *V)
820: {
824: if (ksp->pc_side == PC_RIGHT) {
825: if (ksp->pc) {
826: if (v) {
827: KSP_PCApply(ksp,ksp->vec_sol,v); *V = v;
828: } else SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Not working with right preconditioner");
829: } else {
830: if (v) {
831: VecCopy(ksp->vec_sol,v); *V = v;
832: } else *V = ksp->vec_sol;
833: }
834: } else if (ksp->pc_side == PC_SYMMETRIC) {
835: if (ksp->pc) {
836: if (ksp->transpose_solve) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Not working with symmetric preconditioner and transpose solve");
837: if (v) {
838: PCApplySymmetricRight(ksp->pc,ksp->vec_sol,v);
839: *V = v;
840: } else SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Not working with symmetric preconditioner");
841: } else {
842: if (v) {
843: VecCopy(ksp->vec_sol,v); *V = v;
844: } else *V = ksp->vec_sol;
845: }
846: } else {
847: if (v) {
848: VecCopy(ksp->vec_sol,v); *V = v;
849: } else *V = ksp->vec_sol;
850: }
851: return(0);
852: }
856: /*
857: KSPBuildResidualDefault - Default code to compute the residual.
859: Input Parameters:
860: . ksp - iterative context
861: . t - pointer to temporary vector
862: . v - pointer to user vector
864: Output Parameter:
865: . V - pointer to a vector containing the residual
867: Level: advanced
869: Developers Note: This is PETSC_EXTERN because it may be used by user written plugin KSP implementations
871: .keywords: KSP, build, residual, default
873: .seealso: KSPBuildSolutionDefault()
874: */
875: PetscErrorCode KSPBuildResidualDefault(KSP ksp,Vec t,Vec v,Vec *V)
876: {
878: MatStructure pflag;
879: Mat Amat,Pmat;
882: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
883: PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);
884: KSPBuildSolution(ksp,t,NULL);
885: KSP_MatMult(ksp,Amat,t,v);
886: VecAYPX(v,-1.0,ksp->vec_rhs);
887: *V = v;
888: return(0);
889: }
893: /*@C
894: KSPGetVecs - Gets a number of work vectors.
896: Input Parameters:
897: + ksp - iterative context
898: . rightn - number of right work vectors
899: - leftn - number of left work vectors to allocate
901: Output Parameter:
902: + right - the array of vectors created
903: - left - the array of left vectors
905: Note: The right vector has as many elements as the matrix has columns. The left
906: vector has as many elements as the matrix has rows.
908: Level: advanced
910: .seealso: MatGetVecs()
912: @*/
913: PetscErrorCode KSPGetVecs(KSP ksp,PetscInt rightn, Vec **right,PetscInt leftn,Vec **left)
914: {
916: Vec vecr,vecl;
919: if (rightn) {
920: if (!right) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_INCOMP,"You asked for right vectors but did not pass a pointer to hold them");
921: if (ksp->vec_sol) vecr = ksp->vec_sol;
922: else {
923: if (ksp->dm) {
924: DMGetGlobalVector(ksp->dm,&vecr);
925: } else {
926: Mat mat;
927: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
928: PCGetOperators(ksp->pc,&mat,NULL,NULL);
929: MatGetVecs(mat,&vecr,NULL);
930: }
931: }
932: VecDuplicateVecs(vecr,rightn,right);
933: if (!ksp->vec_sol) {
934: if (ksp->dm) {
935: DMRestoreGlobalVector(ksp->dm,&vecr);
936: } else {
937: VecDestroy(&vecr);
938: }
939: }
940: }
941: if (leftn) {
942: if (!left) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_INCOMP,"You asked for left vectors but did not pass a pointer to hold them");
943: if (ksp->vec_rhs) vecl = ksp->vec_rhs;
944: else {
945: if (ksp->dm) {
946: DMGetGlobalVector(ksp->dm,&vecl);
947: } else {
948: Mat mat;
949: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
950: PCGetOperators(ksp->pc,&mat,NULL,NULL);
951: MatGetVecs(mat,NULL,&vecl);
952: }
953: }
954: VecDuplicateVecs(vecl,leftn,left);
955: if (!ksp->vec_rhs) {
956: if (ksp->dm) {
957: DMRestoreGlobalVector(ksp->dm,&vecl);
958: } else {
959: VecDestroy(&vecl);
960: }
961: }
962: }
963: return(0);
964: }
968: /*
969: KSPSetWorkVecs - Sets a number of work vectors into a KSP object
971: Input Parameters:
972: . ksp - iterative context
973: . nw - number of work vectors to allocate
975: Developers Note: This is PETSC_EXTERN because it may be used by user written plugin KSP implementations
977: */
978: PetscErrorCode KSPSetWorkVecs(KSP ksp,PetscInt nw)
979: {
983: VecDestroyVecs(ksp->nwork,&ksp->work);
984: ksp->nwork = nw;
985: KSPGetVecs(ksp,nw,&ksp->work,0,NULL);
986: PetscLogObjectParents(ksp,nw,ksp->work);
987: return(0);
988: }
992: /*
993: KSPDestroyDefault - Destroys a iterative context variable for methods with
994: no separate context. Preferred calling sequence KSPDestroy().
996: Input Parameter:
997: . ksp - the iterative context
999: Developers Note: This is PETSC_EXTERN because it may be used by user written plugin KSP implementations
1001: */
1002: PetscErrorCode KSPDestroyDefault(KSP ksp)
1003: {
1008: PetscFree(ksp->data);
1009: return(0);
1010: }
1014: /*@
1015: KSPGetConvergedReason - Gets the reason the KSP iteration was stopped.
1017: Not Collective
1019: Input Parameter:
1020: . ksp - the KSP context
1022: Output Parameter:
1023: . reason - negative value indicates diverged, positive value converged, see KSPConvergedReason
1025: Possible values for reason:
1026: + KSP_CONVERGED_RTOL (residual 2-norm decreased by a factor of rtol, from 2-norm of right hand side)
1027: . KSP_CONVERGED_ATOL (residual 2-norm less than abstol)
1028: . KSP_CONVERGED_ITS (used by the preonly preconditioner that always uses ONE iteration, or when the KSPSkipConverged() convergence
1029: test routine is set.
1030: . KSP_CONVERGED_CG_NEG_CURVE
1031: . KSP_CONVERGED_CG_CONSTRAINED
1032: . KSP_CONVERGED_STEP_LENGTH
1033: . KSP_DIVERGED_ITS (required more than its to reach convergence)
1034: . KSP_DIVERGED_DTOL (residual norm increased by a factor of divtol)
1035: . KSP_DIVERGED_NANORINF (residual norm became Not-a-number or Inf likely due to 0/0)
1036: . KSP_DIVERGED_BREAKDOWN (generic breakdown in method)
1037: - KSP_DIVERGED_BREAKDOWN_BICG (Initial residual is orthogonal to preconditioned initial
1038: residual. Try a different preconditioner, or a different initial Level.)
1040: See also manual page for each reason.
1042: guess: beginner
1044: Notes: Can only be called after the call the KSPSolve() is complete.
1046: Level: intermediate
1048: .keywords: KSP, nonlinear, set, convergence, test
1050: .seealso: KSPSetConvergenceTest(), KSPDefaultConverged(), KSPSetTolerances(), KSPConvergedReason
1051: @*/
1052: PetscErrorCode KSPGetConvergedReason(KSP ksp,KSPConvergedReason *reason)
1053: {
1057: *reason = ksp->reason;
1058: return(0);
1059: }
1061: #include <petsc-private/dmimpl.h>
1064: /*@
1065: KSPSetDM - Sets the DM that may be used by some preconditioners
1067: Logically Collective on KSP
1069: Input Parameters:
1070: + ksp - the preconditioner context
1071: - dm - the dm
1073: Notes: If this is used then the KSP will attempt to use the DM to create the matrix and use the routine
1074: set with DMKSPSetComputeOperators(). Use KSPSetDMActive(ksp,PETSC_FALSE) to instead use the matrix
1075: you've provided with KSPSetOperators().
1077: Level: intermediate
1079: .seealso: KSPGetDM(), KSPSetDMActive(), KSPSetComputeOperators(), KSPSetComputeRHS(), KSPSetComputeInitialGuess(), DMKSPSetComputeOperators(), DMKSPSetComputeRHS(), DMKSPSetComputeInitialGuess()
1080: @*/
1081: PetscErrorCode KSPSetDM(KSP ksp,DM dm)
1082: {
1084: PC pc;
1088: if (dm) {PetscObjectReference((PetscObject)dm);}
1089: if (ksp->dm) { /* Move the DMSNES context over to the new DM unless the new DM already has one */
1090: if (ksp->dm->dmksp && ksp->dmAuto && !dm->dmksp) {
1091: DMKSP kdm;
1092: DMCopyDMKSP(ksp->dm,dm);
1093: DMGetDMKSP(ksp->dm,&kdm);
1094: if (kdm->originaldm == ksp->dm) kdm->originaldm = dm; /* Grant write privileges to the replacement DM */
1095: }
1096: DMDestroy(&ksp->dm);
1097: }
1098: ksp->dm = dm;
1099: ksp->dmAuto = PETSC_FALSE;
1100: KSPGetPC(ksp,&pc);
1101: PCSetDM(pc,dm);
1102: ksp->dmActive = PETSC_TRUE;
1103: return(0);
1104: }
1108: /*@
1109: KSPSetDMActive - Indicates the DM should be used to generate the linear system matrix and right hand side
1111: Logically Collective on KSP
1113: Input Parameters:
1114: + ksp - the preconditioner context
1115: - flg - use the DM
1117: Notes:
1118: By default KSPSetDM() sets the DM as active, call KSPSetDMActive(ksp,PETSC_FALSE); after KSPSetDM(ksp,dm) to not have the KSP object use the DM to generate the matrices.
1120: Level: intermediate
1122: .seealso: KSPGetDM(), KSPSetDM(), SNESSetDM(), KSPSetComputeOperators(), KSPSetComputeRHS(), KSPSetComputeInitialGuess()
1123: @*/
1124: PetscErrorCode KSPSetDMActive(KSP ksp,PetscBool flg)
1125: {
1129: ksp->dmActive = flg;
1130: return(0);
1131: }
1135: /*@
1136: KSPGetDM - Gets the DM that may be used by some preconditioners
1138: Not Collective
1140: Input Parameter:
1141: . ksp - the preconditioner context
1143: Output Parameter:
1144: . dm - the dm
1146: Level: intermediate
1149: .seealso: KSPSetDM(), KSPSetDMActive()
1150: @*/
1151: PetscErrorCode KSPGetDM(KSP ksp,DM *dm)
1152: {
1157: if (!ksp->dm) {
1158: DMShellCreate(PetscObjectComm((PetscObject)ksp),&ksp->dm);
1159: ksp->dmAuto = PETSC_TRUE;
1160: }
1161: *dm = ksp->dm;
1162: return(0);
1163: }
1167: /*@
1168: KSPSetApplicationContext - Sets the optional user-defined context for the linear solver.
1170: Logically Collective on KSP
1172: Input Parameters:
1173: + ksp - the KSP context
1174: - usrP - optional user context
1176: Level: intermediate
1178: .keywords: KSP, set, application, context
1180: .seealso: KSPGetApplicationContext()
1181: @*/
1182: PetscErrorCode KSPSetApplicationContext(KSP ksp,void *usrP)
1183: {
1185: PC pc;
1189: ksp->user = usrP;
1190: KSPGetPC(ksp,&pc);
1191: PCSetApplicationContext(pc,usrP);
1192: return(0);
1193: }
1197: /*@
1198: KSPGetApplicationContext - Gets the user-defined context for the linear solver.
1200: Not Collective
1202: Input Parameter:
1203: . ksp - KSP context
1205: Output Parameter:
1206: . usrP - user context
1208: Level: intermediate
1210: .keywords: KSP, get, application, context
1212: .seealso: KSPSetApplicationContext()
1213: @*/
1214: PetscErrorCode KSPGetApplicationContext(KSP ksp,void *usrP)
1215: {
1218: *(void**)usrP = ksp->user;
1219: return(0);
1220: }