Actual source code: lsqr.c
1: /* lourens.vanzanen@shell.com contributed the standard error estimates of the solution, Jul 25, 2006 */
2: /* Bas van't Hof contributed the preconditioned aspects Feb 10, 2010 */
4: #define SWAP(a, b, c) \
5: do { \
6: c = a; \
7: a = b; \
8: b = c; \
9: } while (0)
11: #include <petsc/private/kspimpl.h>
12: #include <petscdraw.h>
14: typedef struct {
15: PetscInt nwork_n, nwork_m;
16: Vec *vwork_m; /* work vectors of length m, where the system is size m x n */
17: Vec *vwork_n; /* work vectors of length n */
18: Vec se; /* Optional standard error vector */
19: PetscBool se_flg; /* flag for -ksp_lsqr_set_standard_error */
20: PetscBool exact_norm; /* flag for -ksp_lsqr_exact_mat_norm */
21: PetscReal arnorm; /* Good estimate of norm((A*inv(Pmat))'*r), where r = A*x - b, used in specific stopping criterion */
22: PetscReal anorm; /* Poor estimate of norm(A*inv(Pmat),'fro') used in specific stopping criterion */
23: /* Backup previous convergence test */
24: PetscErrorCode (*converged)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
25: PetscErrorCode (*convergeddestroy)(void *);
26: void *cnvP;
27: } KSP_LSQR;
29: static PetscErrorCode VecSquare(Vec v)
30: {
31: PetscScalar *x;
32: PetscInt i, n;
34: PetscFunctionBegin;
35: PetscCall(VecGetLocalSize(v, &n));
36: PetscCall(VecGetArray(v, &x));
37: for (i = 0; i < n; i++) x[i] *= PetscConj(x[i]);
38: PetscCall(VecRestoreArray(v, &x));
39: PetscFunctionReturn(PETSC_SUCCESS);
40: }
42: static PetscErrorCode KSPSetUp_LSQR(KSP ksp)
43: {
44: KSP_LSQR *lsqr = (KSP_LSQR *)ksp->data;
45: PetscBool nopreconditioner;
47: PetscFunctionBegin;
48: PetscCall(PetscObjectTypeCompare((PetscObject)ksp->pc, PCNONE, &nopreconditioner));
50: if (lsqr->vwork_m) PetscCall(VecDestroyVecs(lsqr->nwork_m, &lsqr->vwork_m));
52: if (lsqr->vwork_n) PetscCall(VecDestroyVecs(lsqr->nwork_n, &lsqr->vwork_n));
54: lsqr->nwork_m = 2;
55: if (nopreconditioner) lsqr->nwork_n = 4;
56: else lsqr->nwork_n = 5;
57: PetscCall(KSPCreateVecs(ksp, lsqr->nwork_n, &lsqr->vwork_n, lsqr->nwork_m, &lsqr->vwork_m));
59: if (lsqr->se_flg && !lsqr->se) {
60: PetscCall(VecDuplicate(lsqr->vwork_n[0], &lsqr->se));
61: PetscCall(VecSet(lsqr->se, PETSC_INFINITY));
62: } else if (!lsqr->se_flg) {
63: PetscCall(VecDestroy(&lsqr->se));
64: }
65: PetscFunctionReturn(PETSC_SUCCESS);
66: }
68: static PetscErrorCode KSPSolve_LSQR(KSP ksp)
69: {
70: PetscInt i, size1, size2;
71: PetscScalar rho, rhobar, phi, phibar, theta, c, s, tmp, tau;
72: PetscReal beta, alpha, rnorm;
73: Vec X, B, V, V1, U, U1, TMP, W, W2, Z = NULL;
74: Mat Amat, Pmat;
75: KSP_LSQR *lsqr = (KSP_LSQR *)ksp->data;
76: PetscBool diagonalscale, nopreconditioner;
78: PetscFunctionBegin;
79: PetscCall(PCGetDiagonalScale(ksp->pc, &diagonalscale));
80: PetscCheck(!diagonalscale, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Krylov method %s does not support diagonal scaling", ((PetscObject)ksp)->type_name);
82: PetscCall(PCGetOperators(ksp->pc, &Amat, &Pmat));
83: PetscCall(PetscObjectTypeCompare((PetscObject)ksp->pc, PCNONE, &nopreconditioner));
85: /* vectors of length m, where system size is mxn */
86: B = ksp->vec_rhs;
87: U = lsqr->vwork_m[0];
88: U1 = lsqr->vwork_m[1];
90: /* vectors of length n */
91: X = ksp->vec_sol;
92: W = lsqr->vwork_n[0];
93: V = lsqr->vwork_n[1];
94: V1 = lsqr->vwork_n[2];
95: W2 = lsqr->vwork_n[3];
96: if (!nopreconditioner) Z = lsqr->vwork_n[4];
98: /* standard error vector */
99: if (lsqr->se) PetscCall(VecSet(lsqr->se, 0.0));
101: /* Compute initial residual, temporarily use work vector u */
102: if (!ksp->guess_zero) {
103: PetscCall(KSP_MatMult(ksp, Amat, X, U)); /* u <- b - Ax */
104: PetscCall(VecAYPX(U, -1.0, B));
105: } else {
106: PetscCall(VecCopy(B, U)); /* u <- b (x is 0) */
107: }
109: /* Test for nothing to do */
110: PetscCall(VecNorm(U, NORM_2, &rnorm));
111: KSPCheckNorm(ksp, rnorm);
112: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
113: ksp->its = 0;
114: ksp->rnorm = rnorm;
115: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
116: PetscCall(KSPLogResidualHistory(ksp, rnorm));
117: PetscCall(KSPMonitor(ksp, 0, rnorm));
118: PetscCall((*ksp->converged)(ksp, 0, rnorm, &ksp->reason, ksp->cnvP));
119: if (ksp->reason) PetscFunctionReturn(PETSC_SUCCESS);
121: beta = rnorm;
122: PetscCall(VecScale(U, 1.0 / beta));
123: PetscCall(KSP_MatMultHermitianTranspose(ksp, Amat, U, V));
124: if (nopreconditioner) {
125: PetscCall(VecNorm(V, NORM_2, &alpha));
126: KSPCheckNorm(ksp, rnorm);
127: } else {
128: /* this is an application of the preconditioner for the normal equations; not the operator, see the manual page */
129: PetscCall(PCApply(ksp->pc, V, Z));
130: PetscCall(VecDotRealPart(V, Z, &alpha));
131: if (alpha <= 0.0) {
132: ksp->reason = KSP_DIVERGED_BREAKDOWN;
133: PetscFunctionReturn(PETSC_SUCCESS);
134: }
135: alpha = PetscSqrtReal(alpha);
136: PetscCall(VecScale(Z, 1.0 / alpha));
137: }
138: PetscCall(VecScale(V, 1.0 / alpha));
140: if (nopreconditioner) {
141: PetscCall(VecCopy(V, W));
142: } else {
143: PetscCall(VecCopy(Z, W));
144: }
146: if (lsqr->exact_norm) {
147: PetscCall(MatNorm(Amat, NORM_FROBENIUS, &lsqr->anorm));
148: } else lsqr->anorm = 0.0;
150: lsqr->arnorm = alpha * beta;
151: phibar = beta;
152: rhobar = alpha;
153: i = 0;
154: do {
155: if (nopreconditioner) {
156: PetscCall(KSP_MatMult(ksp, Amat, V, U1));
157: } else {
158: PetscCall(KSP_MatMult(ksp, Amat, Z, U1));
159: }
160: PetscCall(VecAXPY(U1, -alpha, U));
161: PetscCall(VecNorm(U1, NORM_2, &beta));
162: KSPCheckNorm(ksp, beta);
163: if (beta > 0.0) {
164: PetscCall(VecScale(U1, 1.0 / beta)); /* beta*U1 = Amat*V - alpha*U */
165: if (!lsqr->exact_norm) lsqr->anorm = PetscSqrtReal(PetscSqr(lsqr->anorm) + PetscSqr(alpha) + PetscSqr(beta));
166: }
168: PetscCall(KSP_MatMultHermitianTranspose(ksp, Amat, U1, V1));
169: PetscCall(VecAXPY(V1, -beta, V));
170: if (nopreconditioner) {
171: PetscCall(VecNorm(V1, NORM_2, &alpha));
172: KSPCheckNorm(ksp, alpha);
173: } else {
174: PetscCall(PCApply(ksp->pc, V1, Z));
175: PetscCall(VecDotRealPart(V1, Z, &alpha));
176: if (alpha <= 0.0) {
177: ksp->reason = KSP_DIVERGED_BREAKDOWN;
178: break;
179: }
180: alpha = PetscSqrtReal(alpha);
181: PetscCall(VecScale(Z, 1.0 / alpha));
182: }
183: PetscCall(VecScale(V1, 1.0 / alpha)); /* alpha*V1 = Amat^T*U1 - beta*V */
184: rho = PetscSqrtScalar(rhobar * rhobar + beta * beta);
185: c = rhobar / rho;
186: s = beta / rho;
187: theta = s * alpha;
188: rhobar = -c * alpha;
189: phi = c * phibar;
190: phibar = s * phibar;
191: tau = s * phi;
193: PetscCall(VecAXPY(X, phi / rho, W)); /* x <- x + (phi/rho) w */
195: if (lsqr->se) {
196: PetscCall(VecCopy(W, W2));
197: PetscCall(VecSquare(W2));
198: PetscCall(VecScale(W2, 1.0 / (rho * rho)));
199: PetscCall(VecAXPY(lsqr->se, 1.0, W2)); /* lsqr->se <- lsqr->se + (w^2/rho^2) */
200: }
201: if (nopreconditioner) {
202: PetscCall(VecAYPX(W, -theta / rho, V1)); /* w <- v - (theta/rho) w */
203: } else {
204: PetscCall(VecAYPX(W, -theta / rho, Z)); /* w <- z - (theta/rho) w */
205: }
207: lsqr->arnorm = alpha * PetscAbsScalar(tau);
208: rnorm = PetscRealPart(phibar);
210: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
211: ksp->its++;
212: ksp->rnorm = rnorm;
213: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
214: PetscCall(KSPLogResidualHistory(ksp, rnorm));
215: PetscCall(KSPMonitor(ksp, i + 1, rnorm));
216: PetscCall((*ksp->converged)(ksp, i + 1, rnorm, &ksp->reason, ksp->cnvP));
217: if (ksp->reason) break;
218: SWAP(U1, U, TMP);
219: SWAP(V1, V, TMP);
221: i++;
222: } while (i < ksp->max_it);
223: if (i >= ksp->max_it && !ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
225: /* Finish off the standard error estimates */
226: if (lsqr->se) {
227: tmp = 1.0;
228: PetscCall(MatGetSize(Amat, &size1, &size2));
229: if (size1 > size2) tmp = size1 - size2;
230: tmp = rnorm / PetscSqrtScalar(tmp);
231: PetscCall(VecSqrtAbs(lsqr->se));
232: PetscCall(VecScale(lsqr->se, tmp));
233: }
234: PetscFunctionReturn(PETSC_SUCCESS);
235: }
237: static PetscErrorCode KSPDestroy_LSQR(KSP ksp)
238: {
239: KSP_LSQR *lsqr = (KSP_LSQR *)ksp->data;
241: PetscFunctionBegin;
242: /* Free work vectors */
243: if (lsqr->vwork_n) PetscCall(VecDestroyVecs(lsqr->nwork_n, &lsqr->vwork_n));
244: if (lsqr->vwork_m) PetscCall(VecDestroyVecs(lsqr->nwork_m, &lsqr->vwork_m));
245: PetscCall(VecDestroy(&lsqr->se));
246: /* Revert convergence test */
247: PetscCall(KSPSetConvergenceTest(ksp, lsqr->converged, lsqr->cnvP, lsqr->convergeddestroy));
248: /* Free the KSP_LSQR context */
249: PetscCall(PetscFree(ksp->data));
250: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPLSQRMonitorResidual_C", NULL));
251: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPLSQRMonitorResidualDrawLG_C", NULL));
252: PetscFunctionReturn(PETSC_SUCCESS);
253: }
255: /*@
256: KSPLSQRSetComputeStandardErrorVec - Compute a vector of standard error estimates during `KSPSolve()` for `KSPLSQR`.
258: Logically Collective
260: Input Parameters:
261: + ksp - iterative context
262: - flg - compute the vector of standard estimates or not
264: Level: intermediate
266: Developer Notes:
267: Vaclav: I'm not sure whether this vector is useful for anything.
269: .seealso: [](ch_ksp), `KSPSolve()`, `KSPLSQR`, `KSPLSQRGetStandardErrorVec()`
270: @*/
271: PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP ksp, PetscBool flg)
272: {
273: KSP_LSQR *lsqr = (KSP_LSQR *)ksp->data;
275: PetscFunctionBegin;
276: lsqr->se_flg = flg;
277: PetscFunctionReturn(PETSC_SUCCESS);
278: }
280: /*@
281: KSPLSQRSetExactMatNorm - Compute exact matrix norm instead of iteratively refined estimate.
283: Not Collective
285: Input Parameters:
286: + ksp - iterative context
287: - flg - compute exact matrix norm or not
289: Level: intermediate
291: Notes:
292: By default, `flg` = `PETSC_FALSE`. This is usually preferred to avoid possibly expensive computation of the norm.
293: For `flg` = `PETSC_TRUE`, we call `MatNorm`(Amat,`NORM_FROBENIUS`,&lsqr->anorm) which will work only for some types of explicitly assembled matrices.
294: This can affect convergence rate as `KSPLSQRConvergedDefault()` assumes different value of ||A|| used in normal equation stopping criterion.
296: .seealso: [](ch_ksp), `KSPSolve()`, `KSPLSQR`, `KSPLSQRGetNorms()`, `KSPLSQRConvergedDefault()`
297: @*/
298: PetscErrorCode KSPLSQRSetExactMatNorm(KSP ksp, PetscBool flg)
299: {
300: KSP_LSQR *lsqr = (KSP_LSQR *)ksp->data;
302: PetscFunctionBegin;
303: lsqr->exact_norm = flg;
304: PetscFunctionReturn(PETSC_SUCCESS);
305: }
307: /*@
308: KSPLSQRGetStandardErrorVec - Get vector of standard error estimates.
309: Only available if -ksp_lsqr_set_standard_error was set to true
310: or `KSPLSQRSetComputeStandardErrorVec`(ksp, `PETSC_TRUE`) was called.
311: Otherwise returns `NULL`.
313: Not Collective
315: Input Parameter:
316: . ksp - iterative context
318: Output Parameter:
319: . se - vector of standard estimates
321: Level: intermediate
323: Developer Notes:
324: Vaclav: I'm not sure whether this vector is useful for anything.
326: .seealso: [](ch_ksp), `KSPSolve()`, `KSPLSQR`, `KSPLSQRSetComputeStandardErrorVec()`
327: @*/
328: PetscErrorCode KSPLSQRGetStandardErrorVec(KSP ksp, Vec *se)
329: {
330: KSP_LSQR *lsqr = (KSP_LSQR *)ksp->data;
332: PetscFunctionBegin;
333: *se = lsqr->se;
334: PetscFunctionReturn(PETSC_SUCCESS);
335: }
337: /*@
338: KSPLSQRGetNorms - Get the norm estimates that `KSPLSQR` computes internally during `KSPSolve()`.
340: Not Collective
342: Input Parameter:
343: . ksp - iterative context
345: Output Parameters:
346: + arnorm - good estimate of $\|(A*Pmat^{-T})*r\|$, where $r = A*x - b$, used in specific stopping criterion
347: - anorm - poor estimate of $\|A*Pmat^{-T}\|_{frobenius}$ used in specific stopping criterion
349: Level: intermediate
351: Notes:
352: Output parameters are meaningful only after `KSPSolve()`.
354: These are the same quantities as normar and norma in MATLAB's `lsqr()`, whose output lsvec is a vector of normar / norma for all iterations.
356: If -ksp_lsqr_exact_mat_norm is set or `KSPLSQRSetExactMatNorm`(ksp, `PETSC_TRUE`) called, then anorm is the exact Frobenius norm.
358: .seealso: [](ch_ksp), `KSPSolve()`, `KSPLSQR`, `KSPLSQRSetExactMatNorm()`
359: @*/
360: PetscErrorCode KSPLSQRGetNorms(KSP ksp, PetscReal *arnorm, PetscReal *anorm)
361: {
362: KSP_LSQR *lsqr = (KSP_LSQR *)ksp->data;
364: PetscFunctionBegin;
365: if (arnorm) *arnorm = lsqr->arnorm;
366: if (anorm) *anorm = lsqr->anorm;
367: PetscFunctionReturn(PETSC_SUCCESS);
368: }
370: static PetscErrorCode KSPLSQRMonitorResidual_LSQR(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
371: {
372: KSP_LSQR *lsqr = (KSP_LSQR *)ksp->data;
373: PetscViewer viewer = vf->viewer;
374: PetscViewerFormat format = vf->format;
375: char normtype[256];
376: PetscInt tablevel;
377: const char *prefix;
379: PetscFunctionBegin;
380: PetscCall(PetscObjectGetTabLevel((PetscObject)ksp, &tablevel));
381: PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ksp, &prefix));
382: PetscCall(PetscStrncpy(normtype, KSPNormTypes[ksp->normtype], sizeof(normtype)));
383: PetscCall(PetscStrtolower(normtype));
384: PetscCall(PetscViewerPushFormat(viewer, format));
385: PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
386: if (n == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, " Residual norm, norm of normal equations, and matrix norm for %s solve.\n", prefix));
387: if (!n) {
388: PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " KSP resid norm %14.12e\n", n, (double)rnorm));
389: } else {
390: PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " KSP resid norm %14.12e normal eq resid norm %14.12e matrix norm %14.12e\n", n, (double)rnorm, (double)lsqr->arnorm, (double)lsqr->anorm));
391: }
392: PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
393: PetscCall(PetscViewerPopFormat(viewer));
394: PetscFunctionReturn(PETSC_SUCCESS);
395: }
397: /*@C
398: KSPLSQRMonitorResidual - Prints the residual norm, as well as the normal equation residual norm, at each iteration of an iterative solver for the `KSPLSQR` solver
400: Collective
402: Input Parameters:
403: + ksp - iterative context
404: . n - iteration number
405: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
406: - vf - The viewer context
408: Options Database Key:
409: . -ksp_lsqr_monitor - Activates `KSPLSQRMonitorResidual()`
411: Level: intermediate
413: .seealso: [](ch_ksp), `KSPLSQR`, `KSPMonitorSet()`, `KSPMonitorResidual()`, `KSPMonitorTrueResidualMaxNorm()`, `KSPLSQRMonitorResidualDrawLG()`
414: @*/
415: PetscErrorCode KSPLSQRMonitorResidual(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
416: {
417: PetscFunctionBegin;
419: PetscAssertPointer(vf, 4);
421: PetscTryMethod(ksp, "KSPLSQRMonitorResidual_C", (KSP, PetscInt, PetscReal, PetscViewerAndFormat *), (ksp, n, rnorm, vf));
422: PetscFunctionReturn(PETSC_SUCCESS);
423: }
425: static PetscErrorCode KSPLSQRMonitorResidualDrawLG_LSQR(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
426: {
427: KSP_LSQR *lsqr = (KSP_LSQR *)ksp->data;
428: PetscViewer viewer = vf->viewer;
429: PetscViewerFormat format = vf->format;
430: PetscDrawLG lg = vf->lg;
431: KSPConvergedReason reason;
432: PetscReal x[2], y[2];
434: PetscFunctionBegin;
435: PetscCall(PetscViewerPushFormat(viewer, format));
436: if (!n) PetscCall(PetscDrawLGReset(lg));
437: x[0] = (PetscReal)n;
438: if (rnorm > 0.0) y[0] = PetscLog10Real(rnorm);
439: else y[0] = -15.0;
440: x[1] = (PetscReal)n;
441: if (lsqr->arnorm > 0.0) y[1] = PetscLog10Real(lsqr->arnorm);
442: else y[1] = -15.0;
443: PetscCall(PetscDrawLGAddPoint(lg, x, y));
444: PetscCall(KSPGetConvergedReason(ksp, &reason));
445: if (n <= 20 || !(n % 5) || reason) {
446: PetscCall(PetscDrawLGDraw(lg));
447: PetscCall(PetscDrawLGSave(lg));
448: }
449: PetscCall(PetscViewerPopFormat(viewer));
450: PetscFunctionReturn(PETSC_SUCCESS);
451: }
453: /*@C
454: KSPLSQRMonitorResidualDrawLG - Plots the true residual norm at each iteration of an iterative solver for the `KSPLSQR` solver
456: Collective
458: Input Parameters:
459: + ksp - iterative context
460: . n - iteration number
461: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
462: - vf - The viewer context
464: Options Database Key:
465: . -ksp_lsqr_monitor draw::draw_lg - Activates `KSPMonitorTrueResidualDrawLG()`
467: Level: intermediate
469: .seealso: [](ch_ksp), `KSPLSQR`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `KSPLSQRMonitorResidual()`, `KSPLSQRMonitorResidualDrawLGCreate()`
470: @*/
471: PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
472: {
473: PetscFunctionBegin;
475: PetscAssertPointer(vf, 4);
478: PetscTryMethod(ksp, "KSPLSQRMonitorResidualDrawLG_C", (KSP, PetscInt, PetscReal, PetscViewerAndFormat *), (ksp, n, rnorm, vf));
479: PetscFunctionReturn(PETSC_SUCCESS);
480: }
482: /*@C
483: KSPLSQRMonitorResidualDrawLGCreate - Creates the line graph object for the `KSPLSQR` residual and normal equation residual norm
485: Collective
487: Input Parameters:
488: + viewer - The `PetscViewer`
489: . format - The viewer format
490: - ctx - An optional user context
492: Output Parameter:
493: . vf - The `PetscViewerAndFormat`
495: Level: intermediate
497: .seealso: [](ch_ksp), `KSPLSQR`, `KSPMonitorSet()`, `KSPLSQRMonitorResidual()`, `KSPLSQRMonitorResidualDrawLG()`
498: @*/
499: PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer viewer, PetscViewerFormat format, void *ctx, PetscViewerAndFormat **vf)
500: {
501: const char *names[] = {"residual", "normal eqn residual"};
503: PetscFunctionBegin;
504: PetscCall(PetscViewerAndFormatCreate(viewer, format, vf));
505: (*vf)->data = ctx;
506: PetscCall(KSPMonitorLGCreate(PetscObjectComm((PetscObject)viewer), NULL, NULL, "Log Residual Norm", 2, names, PETSC_DECIDE, PETSC_DECIDE, 400, 300, &(*vf)->lg));
507: PetscFunctionReturn(PETSC_SUCCESS);
508: }
510: static PetscErrorCode KSPSetFromOptions_LSQR(KSP ksp, PetscOptionItems *PetscOptionsObject)
511: {
512: KSP_LSQR *lsqr = (KSP_LSQR *)ksp->data;
514: PetscFunctionBegin;
515: PetscOptionsHeadBegin(PetscOptionsObject, "KSP LSQR Options");
516: PetscCall(PetscOptionsBool("-ksp_lsqr_compute_standard_error", "Set Standard Error Estimates of Solution", "KSPLSQRSetComputeStandardErrorVec", lsqr->se_flg, &lsqr->se_flg, NULL));
517: PetscCall(PetscOptionsBool("-ksp_lsqr_exact_mat_norm", "Compute exact matrix norm instead of iteratively refined estimate", "KSPLSQRSetExactMatNorm", lsqr->exact_norm, &lsqr->exact_norm, NULL));
518: PetscCall(KSPMonitorSetFromOptions(ksp, "-ksp_lsqr_monitor", "lsqr_residual", NULL));
519: PetscOptionsHeadEnd();
520: PetscFunctionReturn(PETSC_SUCCESS);
521: }
523: static PetscErrorCode KSPView_LSQR(KSP ksp, PetscViewer viewer)
524: {
525: KSP_LSQR *lsqr = (KSP_LSQR *)ksp->data;
526: PetscBool iascii;
528: PetscFunctionBegin;
529: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
530: if (iascii) {
531: if (lsqr->se) {
532: PetscReal rnorm;
533: PetscCall(VecNorm(lsqr->se, NORM_2, &rnorm));
534: PetscCall(PetscViewerASCIIPrintf(viewer, " norm of standard error %g, iterations %" PetscInt_FMT "\n", (double)rnorm, ksp->its));
535: } else {
536: PetscCall(PetscViewerASCIIPrintf(viewer, " standard error not computed\n"));
537: }
538: if (lsqr->exact_norm) {
539: PetscCall(PetscViewerASCIIPrintf(viewer, " using exact matrix norm\n"));
540: } else {
541: PetscCall(PetscViewerASCIIPrintf(viewer, " using inexact matrix norm\n"));
542: }
543: }
544: PetscFunctionReturn(PETSC_SUCCESS);
545: }
547: /*@C
548: KSPLSQRConvergedDefault - Determines convergence of the `KSPLSQR` Krylov method.
550: Collective
552: Input Parameters:
553: + ksp - iterative context
554: . n - iteration number
555: . rnorm - 2-norm residual value (may be estimated)
556: - ctx - convergence context which must have been created by `KSPConvergedDefaultCreate()`
558: Output Parameter:
559: . reason - the convergence reason
561: Level: advanced
563: Notes:
564: This is not called directly but rather is passed to `KSPSetConvergenceTest()`. It is used automatically by `KSPLSQR`
566: `KSPConvergedDefault()` is called first to check for convergence in $A*x=b$.
567: If that does not determine convergence then checks convergence for the least squares problem, i.e. in min{|b-A*x|}.
568: Possible convergence for the least squares problem (which is based on the residual of the normal equations) are `KSP_CONVERGED_RTOL_NORMAL` norm
569: and `KSP_CONVERGED_ATOL_NORMAL`.
571: `KSP_CONVERGED_RTOL_NORMAL` is returned if $||A^T*r|| < rtol * ||A|| * ||r||$.
572: Matrix norm $||A||$ is iteratively refined estimate, see `KSPLSQRGetNorms()`.
573: This criterion is largely compatible with that in MATLAB `lsqr()`.
575: .seealso: [](ch_ksp), `KSPLSQR`, `KSPSetConvergenceTest()`, `KSPSetTolerances()`, `KSPConvergedSkip()`, `KSPConvergedReason`, `KSPGetConvergedReason()`,
576: `KSPConvergedDefaultSetUIRNorm()`, `KSPConvergedDefaultSetUMIRNorm()`, `KSPConvergedDefaultCreate()`, `KSPConvergedDefaultDestroy()`, `KSPConvergedDefault()`, `KSPLSQRGetNorms()`, `KSPLSQRSetExactMatNorm()`
577: @*/
578: PetscErrorCode KSPLSQRConvergedDefault(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, void *ctx)
579: {
580: KSP_LSQR *lsqr = (KSP_LSQR *)ksp->data;
582: PetscFunctionBegin;
583: /* check for convergence in A*x=b */
584: PetscCall(KSPConvergedDefault(ksp, n, rnorm, reason, ctx));
585: if (!n || *reason) PetscFunctionReturn(PETSC_SUCCESS);
587: /* check for convergence in min{|b-A*x|} */
588: if (lsqr->arnorm < ksp->abstol) {
589: PetscCall(PetscInfo(ksp, "LSQR solver has converged. Normal equation residual %14.12e is less than absolute tolerance %14.12e at iteration %" PetscInt_FMT "\n", (double)lsqr->arnorm, (double)ksp->abstol, n));
590: *reason = KSP_CONVERGED_ATOL_NORMAL;
591: } else if (lsqr->arnorm < ksp->rtol * lsqr->anorm * rnorm) {
592: PetscCall(PetscInfo(ksp, "LSQR solver has converged. Normal equation residual %14.12e is less than rel. tol. %14.12e times %s Frobenius norm of matrix %14.12e times residual %14.12e at iteration %" PetscInt_FMT "\n", (double)lsqr->arnorm,
593: (double)ksp->rtol, lsqr->exact_norm ? "exact" : "approx.", (double)lsqr->anorm, (double)rnorm, n));
594: *reason = KSP_CONVERGED_RTOL_NORMAL;
595: }
596: PetscFunctionReturn(PETSC_SUCCESS);
597: }
599: /*MC
600: KSPLSQR - Implements LSQR {cite}`paige.saunders:lsqr`
602: Options Database Keys:
603: + -ksp_lsqr_set_standard_error - set standard error estimates of solution, see `KSPLSQRSetComputeStandardErrorVec()` and `KSPLSQRGetStandardErrorVec()`
604: . -ksp_lsqr_exact_mat_norm - compute exact matrix norm instead of iteratively refined estimate, see `KSPLSQRSetExactMatNorm()`
605: - -ksp_lsqr_monitor - monitor residual norm, norm of residual of normal equations A'*A x = A' b, and estimate of matrix norm ||A||
607: Level: beginner
609: Notes:
610: Supports non-square (rectangular) matrices.
612: This variant, when applied with no preconditioning is identical to the original algorithm in exact arithmetic; however, in practice, with no preconditioning
613: due to inexact arithmetic, it can converge differently. Hence when no preconditioner is used (`PCType` `PCNONE`) it automatically reverts to the original algorithm.
615: With the PETSc built-in preconditioners, such as `PCICC`, one should call `KSPSetOperators`(ksp,A,A'*A)) since the preconditioner needs to work
616: for the normal equations A'*A.
618: Supports only left preconditioning.
620: For least squares problems with nonzero residual $A*x - b$, there are additional convergence tests for the residual of the normal equations, $A^T*(b - Ax)$, see `KSPLSQRConvergedDefault()`.
622: In exact arithmetic the LSQR method (with no preconditioning) is identical to the `KSPCG` algorithm applied to the normal equations.
623: The preconditioned variant was implemented by Bas van't Hof and is essentially a left preconditioning for the Normal Equations.
624: It appears the implementation with preconditioning tracks the true norm of the residual and uses that in the convergence test.
626: Developer Note:
627: How is this related to the `KSPCGNE` implementation? One difference is that `KSPCGNE` applies
628: the preconditioner transpose times the preconditioner, so one does not need to pass $A^T*A$ as the third argument to `KSPSetOperators()`.
630: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPSolve()`, `KSPLSQRConvergedDefault()`, `KSPLSQRSetComputeStandardErrorVec()`, `KSPLSQRGetStandardErrorVec()`, `KSPLSQRSetExactMatNorm()`, `KSPLSQRMonitorResidualDrawLGCreate()`, `KSPLSQRMonitorResidualDrawLG()`, `KSPLSQRMonitorResidual()`
631: M*/
632: PETSC_EXTERN PetscErrorCode KSPCreate_LSQR(KSP ksp)
633: {
634: KSP_LSQR *lsqr;
635: void *ctx;
637: PetscFunctionBegin;
638: PetscCall(PetscNew(&lsqr));
639: lsqr->se = NULL;
640: lsqr->se_flg = PETSC_FALSE;
641: lsqr->exact_norm = PETSC_FALSE;
642: lsqr->anorm = -1.0;
643: lsqr->arnorm = -1.0;
644: ksp->data = (void *)lsqr;
645: PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_LEFT, 3));
647: ksp->ops->setup = KSPSetUp_LSQR;
648: ksp->ops->solve = KSPSolve_LSQR;
649: ksp->ops->destroy = KSPDestroy_LSQR;
650: ksp->ops->setfromoptions = KSPSetFromOptions_LSQR;
651: ksp->ops->view = KSPView_LSQR;
653: /* Backup current convergence test; remove destroy routine from KSP to prevent destroying the convergence context in KSPSetConvergenceTest() */
654: PetscCall(KSPGetAndClearConvergenceTest(ksp, &lsqr->converged, &lsqr->cnvP, &lsqr->convergeddestroy));
655: /* Override current convergence test */
656: PetscCall(KSPConvergedDefaultCreate(&ctx));
657: PetscCall(KSPSetConvergenceTest(ksp, KSPLSQRConvergedDefault, ctx, KSPConvergedDefaultDestroy));
658: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPLSQRMonitorResidual_C", KSPLSQRMonitorResidual_LSQR));
659: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPLSQRMonitorResidualDrawLG_C", KSPLSQRMonitorResidualDrawLG_LSQR));
660: PetscFunctionReturn(PETSC_SUCCESS);
661: }