Actual source code: iterativ.c
2: /*
3: This file contains some simple default routines.
4: These routines should be SHORT, since they will be included in every
5: executable image that uses the iterative routines (note that, through
6: the registry system, we provide a way to load only the truely necessary
7: files)
8: */
9: #include src/ksp/ksp/kspimpl.h
13: /*
14: KSPDefaultFreeWork - Free work vectors
16: Input Parameters:
17: . ksp - iterative context
18: */
19: PetscErrorCode KSPDefaultFreeWork(KSP ksp)
20: {
24: if (ksp->work) {
25: VecDestroyVecs(ksp->work,ksp->nwork);
26: ksp->work = PETSC_NULL;
27: }
28: return(0);
29: }
33: /*@C
34: KSPGetResidualNorm - Gets the last (approximate preconditioned)
35: residual norm that has been computed.
36:
37: Not Collective
39: Input Parameters:
40: . ksp - the iterative context
42: Output Parameters:
43: . rnorm - residual norm
45: Level: intermediate
47: .keywords: KSP, get, residual norm
49: .seealso: KSPComputeResidual()
50: @*/
51: PetscErrorCode KSPGetResidualNorm(KSP ksp,PetscReal *rnorm)
52: {
56: *rnorm = ksp->rnorm;
57: return(0);
58: }
62: /*@
63: KSPGetIterationNumber - Gets the current iteration number; if the
64: KSPSolve() is complete, returns the number of iterations
65: used.
66:
67: Not Collective
69: Input Parameters:
70: . ksp - the iterative context
72: Output Parameters:
73: . its - number of iterations
75: Level: intermediate
77: Notes:
78: During the ith iteration this returns i-1
79: .keywords: KSP, get, residual norm
81: .seealso: KSPComputeResidual(), KSPGetResidualNorm()
82: @*/
83: PetscErrorCode KSPGetIterationNumber(KSP ksp,PetscInt *its)
84: {
88: *its = ksp->its;
89: return(0);
90: }
94: /*@C
95: KSPSingularValueMonitor - Prints the two norm of the true residual and
96: estimation of the extreme singular values of the preconditioned problem
97: at each iteration.
98:
99: Collective on KSP
101: Input Parameters:
102: + ksp - the iterative context
103: . n - the iteration
104: - rnorm - the two norm of the residual
106: Options Database Key:
107: . -ksp_singmonitor - Activates KSPSingularValueMonitor()
109: Notes:
110: The CG solver uses the Lanczos technique for eigenvalue computation,
111: while GMRES uses the Arnoldi technique; other iterative methods do
112: not currently compute singular values.
114: Level: intermediate
116: .keywords: KSP, CG, default, monitor, extreme, singular values, Lanczos, Arnoldi
118: .seealso: KSPComputeExtremeSingularValues()
119: @*/
120: PetscErrorCode KSPSingularValueMonitor(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy)
121: {
122: PetscReal emin,emax,c;
127: if (!ksp->calc_sings) {
128: PetscPrintf(ksp->comm,"%3D KSP Residual norm %14.12e \n",n,rnorm);
129: } else {
130: KSPComputeExtremeSingularValues(ksp,&emax,&emin);
131: c = emax/emin;
132: PetscPrintf(ksp->comm,"%3D KSP Residual norm %14.12e %% max %g min %g max/min %g\n",n,rnorm,emax,emin,c);
133: }
134: return(0);
135: }
139: /*@C
140: KSPVecViewMonitor - Monitors progress of the KSP solvers by calling
141: VecView() for the approximate solution at each iteration.
143: Collective on KSP
145: Input Parameters:
146: + ksp - the KSP context
147: . its - iteration number
148: . fgnorm - 2-norm of residual (or gradient)
149: - dummy - either a viewer or PETSC_NULL
151: Level: intermediate
153: Notes:
154: For some Krylov methods such as GMRES constructing the solution at
155: each iteration is expensive, hence using this will slow the code.
157: .keywords: KSP, nonlinear, vector, monitor, view
159: .seealso: KSPSetMonitor(), KSPDefaultMonitor(), VecView()
160: @*/
161: PetscErrorCode KSPVecViewMonitor(KSP ksp,PetscInt its,PetscReal fgnorm,void *dummy)
162: {
164: Vec x;
165: PetscViewer viewer = (PetscViewer) dummy;
168: KSPBuildSolution(ksp,PETSC_NULL,&x);
169: if (!viewer) {
170: MPI_Comm comm;
171: PetscObjectGetComm((PetscObject)ksp,&comm);
172: viewer = PETSC_VIEWER_DRAW_(comm);
173: }
174: VecView(x,viewer);
176: return(0);
177: }
181: /*@C
182: KSPDefaultMonitor - Print the residual norm at each iteration of an
183: iterative solver.
185: Collective on KSP
187: Input Parameters:
188: + ksp - iterative context
189: . n - iteration number
190: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
191: - dummy - unused monitor context
193: Level: intermediate
195: .keywords: KSP, default, monitor, residual
197: .seealso: KSPSetMonitor(), KSPTrueMonitor(), KSPLGMonitorCreate()
198: @*/
199: PetscErrorCode KSPDefaultMonitor(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy)
200: {
202: PetscViewer viewer = (PetscViewer) dummy;
205: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(ksp->comm);
206: PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm %14.12e \n",n,rnorm);
207: return(0);
208: }
212: /*@C
213: KSPTrueMonitor - Prints the true residual norm as well as the preconditioned
214: residual norm at each iteration of an iterative solver.
216: Collective on KSP
218: Input Parameters:
219: + ksp - iterative context
220: . n - iteration number
221: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
222: - dummy - unused monitor context
224: Options Database Key:
225: . -ksp_truemonitor - Activates KSPTrueMonitor()
227: Notes:
228: When using right preconditioning, these values are equivalent.
230: When using either ICC or ILU preconditioners in BlockSolve95
231: (via MATMPIROWBS matrix format), then use this monitor will
232: print both the residual norm associated with the original
233: (unscaled) matrix.
235: Level: intermediate
237: .keywords: KSP, default, monitor, residual
239: .seealso: KSPSetMonitor(), KSPDefaultMonitor(), KSPLGMonitorCreate()
240: @*/
241: PetscErrorCode KSPTrueMonitor(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy)
242: {
244: Vec resid,work;
245: PetscReal scnorm;
246: PC pc;
247: Mat A,B;
248: PetscViewer viewer = (PetscViewer) dummy;
249:
251: VecDuplicate(ksp->vec_rhs,&work);
252: KSPBuildResidual(ksp,0,work,&resid);
254: /*
255: Unscale the residual if the matrix is, for example, a BlockSolve matrix
256: but only if both matrices are the same matrix, since only then would
257: they be scaled.
258: */
259: VecCopy(resid,work);
260: KSPGetPC(ksp,&pc);
261: PCGetOperators(pc,&A,&B,PETSC_NULL);
262: if (A == B) {
263: MatUnScaleSystem(A,PETSC_NULL,work);
264: }
265: VecNorm(work,NORM_2,&scnorm);
266: VecDestroy(work);
267: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(ksp->comm);
268: PetscViewerASCIIPrintf(viewer,"%3D KSP preconditioned resid norm %14.12e true resid norm %14.12e\n",n,rnorm,scnorm);
269: return(0);
270: }
274: /*
275: Default (short) KSP Monitor, same as KSPDefaultMonitor() except
276: it prints fewer digits of the residual as the residual gets smaller.
277: This is because the later digits are meaningless and are often
278: different on different machines; by using this routine different
279: machines will usually generate the same output.
280: */
281: PetscErrorCode KSPDefaultSMonitor(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy)
282: {
284: PetscViewer viewer = (PetscViewer) dummy;
287: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(ksp->comm);
289: if (fnorm > 1.e-9) {
290: PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm %g \n",its,fnorm);
291: } else if (fnorm > 1.e-11){
292: PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm %5.3e \n",its,fnorm);
293: } else {
294: PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm < 1.e-11\n",its);
295: }
296: return(0);
297: }
301: /*@C
302: KSPSkipConverged - Convergence test that NEVER returns as converged.
304: Collective on KSP
306: Input Parameters:
307: + ksp - iterative context
308: . n - iteration number
309: . rnorm - 2-norm residual value (may be estimated)
310: - dummy - unused convergence context
312: Returns:
313: . reason - always KSP_CONVERGED_ITERATING
315: Notes:
316: This is used as the convergence test with the option KSPSetNormType(ksp,KSP_NO_NORM),
317: since norms of the residual are not computed. Convergence is then declared
318: after a fixed number of iterations have been used. Useful when one is
319: using CG or Bi-CG-stab as a smoother.
320:
321: Level: advanced
323: .keywords: KSP, default, convergence, residual
325: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPSetNormType()
326: @*/
327: PetscErrorCode KSPSkipConverged(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *dummy)
328: {
331: return(0);
332: }
336: /*@C
337: KSPDefaultConverged - Determines convergence of
338: the iterative solvers (default code).
340: Collective on KSP
342: Input Parameters:
343: + ksp - iterative context
344: . n - iteration number
345: . rnorm - 2-norm residual value (may be estimated)
346: - dummy - unused convergence context
348: Returns:
349: + positive - if the iteration has converged;
350: . negative - if residual norm exceeds divergence threshold;
351: - 0 - otherwise.
353: Notes:
354: KSPDefaultConverged() reaches convergence when
355: $ rnorm < MAX (rtol * rnorm_0, abstol);
356: Divergence is detected if
357: $ rnorm > dtol * rnorm_0,
359: where
360: + rtol = relative tolerance,
361: . abstol = absolute tolerance.
362: . dtol = divergence tolerance,
363: - rnorm_0 = initial residual norm
365: Use KSPSetTolerances() to alter the defaults for rtol, abstol, dtol.
367: The precise values of reason are macros such as KSP_CONVERGED_RTOL, which
368: are defined in petscksp.h.
370: Level: intermediate
372: .keywords: KSP, default, convergence, residual
374: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPSkipConverged(), KSPConvergedReason, KSPGetConvergedReason()
375: @*/
376: PetscErrorCode KSPDefaultConverged(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *dummy)
377: {
381: *reason = KSP_CONVERGED_ITERATING;
383: if (!n) {
384: /* if user gives initial guess need to compute norm of b */
385: if (!ksp->guess_zero) {
386: PetscReal snorm;
388: if (ksp->normtype == KSP_UNPRECONDITIONED_NORM) {
389: VecNorm(ksp->vec_rhs,NORM_2,&snorm); /* <- b'*b */
390: } else {
391: Vec z;
392: VecDuplicate(ksp->vec_rhs,&z);
393: KSP_PCApply(ksp,ksp->vec_rhs,z);
394: if (ksp->normtype == KSP_PRECONDITIONED_NORM) {
395: VecNorm(z,NORM_2,&snorm); /* dp <- b'*B'*B*b */
396: } else if (ksp->normtype == KSP_NATURAL_NORM) {
397: PetscScalar norm;
398: VecDot(ksp->vec_rhs,z,&norm);
399: snorm = sqrt(PetscAbsScalar(norm)); /* dp <- b'*B*b */
400: }
401: VecDestroy(z);
402: }
403: /* handle special case of zero RHS and nonzero guess */
404: if (!snorm) snorm = rnorm;
405: ksp->ttol = PetscMax(ksp->rtol*snorm,ksp->abstol);
406: ksp->rnorm0 = snorm;
407: } else {
408: ksp->ttol = PetscMax(ksp->rtol*rnorm,ksp->abstol);
409: ksp->rnorm0 = rnorm;
410: }
411: }
412: if (rnorm <= ksp->ttol) {
413: if (rnorm < ksp->abstol) {
414: PetscLogInfo(ksp,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g at iteration %D\n",rnorm,ksp->abstol,n);
415: *reason = KSP_CONVERGED_ATOL;
416: } else {
417: PetscLogInfo(ksp,"Linear solver has converged. Residual norm %g is less than relative tolerance %g times initial residual norm %g at iteration %D\n",rnorm,ksp->rtol,ksp->rnorm0,n);
418: *reason = KSP_CONVERGED_RTOL;
419: }
420: } else if (rnorm >= ksp->divtol*ksp->rnorm0) {
421: PetscLogInfo(ksp,"Linear solver is diverging. Initial residual norm %g, current residual norm %g at iteration %D\n",ksp->rnorm0,rnorm,n);
422: *reason = KSP_DIVERGED_DTOL;
423: } else if (rnorm != rnorm) {
424: PetscLogInfo(ksp,"Linear solver has created a not a number (NaN) as the residual norm, declaring divergence\n");
425: *reason = KSP_DIVERGED_DTOL;
426: }
427: return(0);
428: }
432: /*
433: KSPDefaultBuildSolution - Default code to create/move the solution.
435: Input Parameters:
436: + ksp - iterative context
437: - v - pointer to the user's vector
439: Output Parameter:
440: . V - pointer to a vector containing the solution
442: Level: advanced
444: .keywords: KSP, build, solution, default
446: .seealso: KSPGetSolution(), KSPDefaultBuildResidual()
447: */
448: PetscErrorCode KSPDefaultBuildSolution(KSP ksp,Vec v,Vec *V)
449: {
452: if (ksp->pc_side == PC_RIGHT) {
453: if (ksp->pc) {
454: if (v) {KSP_PCApply(ksp,ksp->vec_sol,v); *V = v;}
455: else {SETERRQ(PETSC_ERR_SUP,"Not working with right preconditioner");}
456: } else {
457: if (v) {VecCopy(ksp->vec_sol,v); *V = v;}
458: else { *V = ksp->vec_sol;}
459: }
460: } else if (ksp->pc_side == PC_SYMMETRIC) {
461: if (ksp->pc) {
462: if (ksp->transpose_solve) SETERRQ(PETSC_ERR_SUP,"Not working with symmetric preconditioner and transpose solve");
463: if (v) {PCApplySymmetricRight(ksp->pc,ksp->vec_sol,v); *V = v;}
464: else {SETERRQ(PETSC_ERR_SUP,"Not working with symmetric preconditioner");}
465: } else {
466: if (v) {VecCopy(ksp->vec_sol,v); *V = v;}
467: else { *V = ksp->vec_sol;}
468: }
469: } else {
470: if (v) {VecCopy(ksp->vec_sol,v); *V = v;}
471: else { *V = ksp->vec_sol; }
472: }
473: return(0);
474: }
478: /*
479: KSPDefaultBuildResidual - Default code to compute the residual.
481: Input Parameters:
482: . ksp - iterative context
483: . t - pointer to temporary vector
484: . v - pointer to user vector
486: Output Parameter:
487: . V - pointer to a vector containing the residual
489: Level: advanced
491: .keywords: KSP, build, residual, default
493: .seealso: KSPDefaultBuildSolution()
494: */
495: PetscErrorCode KSPDefaultBuildResidual(KSP ksp,Vec t,Vec v,Vec *V)
496: {
498: MatStructure pflag;
499: Vec T;
500: PetscScalar mone = -1.0;
501: Mat Amat,Pmat;
504: PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);
505: KSPBuildSolution(ksp,t,&T);
506: KSP_MatMult(ksp,Amat,t,v);
507: VecAYPX(&mone,ksp->vec_rhs,v);
508: *V = v;
509: return(0);
510: }
514: /*
515: KSPDefaultGetWork - Gets a number of work vectors.
517: Input Parameters:
518: . ksp - iterative context
519: . nw - number of work vectors to allocate
521: Output Parameter:
522: . work - the array of vectors created
524: */
525: PetscErrorCode KSPGetVecs(KSP ksp,PetscInt nw,Vec **work)
526: {
528: Vec vec;
531: if (ksp->vec_rhs) vec = ksp->vec_rhs;
532: else {
533: Mat pmat;
534: PCGetOperators(ksp->pc,0,&pmat,0);
535: MatGetVecs(pmat,&vec,0);
536: }
537: VecDuplicateVecs(vec,nw,work);
538: if (!ksp->vec_rhs) {
539: VecDestroy(vec);
540: }
541: return(0);
542: }
546: /*
547: KSPDefaultGetWork - Gets a number of work vectors.
549: Input Parameters:
550: . ksp - iterative context
551: . nw - number of work vectors to allocate
553: Notes:
554: Call this only if no work vectors have been allocated
555: */
556: PetscErrorCode KSPDefaultGetWork(KSP ksp,PetscInt nw)
557: {
561: if (ksp->work) {KSPDefaultFreeWork(ksp);}
562: ksp->nwork = nw;
563: KSPGetVecs(ksp,nw,&ksp->work);
564: PetscLogObjectParents(ksp,nw,ksp->work);
565: return(0);
566: }
570: /*
571: KSPDefaultDestroy - Destroys a iterative context variable for methods with
572: no separate context. Preferred calling sequence KSPDestroy().
574: Input Parameter:
575: . ksp - the iterative context
576: */
577: PetscErrorCode KSPDefaultDestroy(KSP ksp)
578: {
583: if (ksp->data) {
584: PetscFree(ksp->data);
585: ksp->data = 0;
586: }
588: /* free work vectors */
589: KSPDefaultFreeWork(ksp);
590: return(0);
591: }
595: /*@C
596: KSPGetConvergedReason - Gets the reason the KSP iteration was stopped.
598: Not Collective
600: Input Parameter:
601: . ksp - the KSP context
603: Output Parameter:
604: . reason - negative value indicates diverged, positive value converged, see KSPConvergedReason
606: Possible values for reason:
607: + KSP_CONVERGED_RTOL (residual 2-norm decreased by a factor of rtol, from 2-norm of right hand side)
608: . KSP_CONVERGED_ATOL (residual 2-norm less than abstol)
609: . KSP_CONVERGED_ITS (used by the preonly preconditioner that always uses ONE iteration)
610: . KSP_CONVERGED_QCG_NEG_CURVE
611: . KSP_CONVERGED_QCG_CONSTRAINED
612: . KSP_CONVERGED_STEP_LENGTH
613: . KSP_DIVERGED_ITS (required more than its to reach convergence)
614: . KSP_DIVERGED_DTOL (residual norm increased by a factor of divtol)
615: . KSP_DIVERGED_BREAKDOWN (generic breakdown in method)
616: - KSP_DIVERGED_BREAKDOWN_BICG (Initial residual is orthogonal to preconditioned initial
617: residual. Try a different preconditioner, or a different initial guess.)
618:
620: Level: beginner
622: Notes: Can only be called after the call the KSPSolve() is complete.
624: .keywords: KSP, nonlinear, set, convergence, test
626: .seealso: KSPSetConvergenceTest(), KSPDefaultConverged(), KSPSetTolerances(), KSPConvergedReason
627: @*/
628: PetscErrorCode KSPGetConvergedReason(KSP ksp,KSPConvergedReason *reason)
629: {
633: *reason = ksp->reason;
634: return(0);
635: }