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:   KSPConvergenceTestFn *converged;
 25:   PetscCtxDestroyFn    *convergeddestroy;
 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:       PetscCall(PetscInfo(ksp, "Diverging due to breakdown alpha (%g) <= 0\n", (double)alpha));
134:       PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "KSPSolve breakdown alpha (%g) <= 0", (double)alpha);
135:       PetscFunctionReturn(PETSC_SUCCESS);
136:     }
137:     alpha = PetscSqrtReal(alpha);
138:     PetscCall(VecScale(Z, 1.0 / alpha));
139:   }
140:   PetscCall(VecScale(V, 1.0 / alpha));

142:   if (nopreconditioner) {
143:     PetscCall(VecCopy(V, W));
144:   } else {
145:     PetscCall(VecCopy(Z, W));
146:   }

148:   if (lsqr->exact_norm) PetscCall(MatNorm(Amat, NORM_FROBENIUS, &lsqr->anorm));
149:   else lsqr->anorm = 0.0;

151:   lsqr->arnorm = alpha * beta;
152:   phibar       = beta;
153:   rhobar       = alpha;
154:   i            = 0;
155:   do {
156:     if (nopreconditioner) {
157:       PetscCall(KSP_MatMult(ksp, Amat, V, U1));
158:     } else {
159:       PetscCall(KSP_MatMult(ksp, Amat, Z, U1));
160:     }
161:     PetscCall(VecAXPY(U1, -alpha, U));
162:     PetscCall(VecNorm(U1, NORM_2, &beta));
163:     KSPCheckNorm(ksp, beta);
164:     if (beta > 0.0) {
165:       PetscCall(VecScale(U1, 1.0 / beta)); /* beta*U1 = Amat*V - alpha*U */
166:       if (!lsqr->exact_norm) lsqr->anorm = PetscSqrtReal(PetscSqr(lsqr->anorm) + PetscSqr(alpha) + PetscSqr(beta));
167:     }

169:     PetscCall(KSP_MatMultHermitianTranspose(ksp, Amat, U1, V1));
170:     PetscCall(VecAXPY(V1, -beta, V));
171:     if (nopreconditioner) {
172:       PetscCall(VecNorm(V1, NORM_2, &alpha));
173:       KSPCheckNorm(ksp, alpha);
174:     } else {
175:       PetscCall(PCApply(ksp->pc, V1, Z));
176:       PetscCall(VecDotRealPart(V1, Z, &alpha));
177:       if (alpha < 0.0) {
178:         PetscCall(PetscInfo(ksp, "Diverging due to breakdown alpha (%g) < 0\n", (double)alpha));
179:         PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "KSPSolve breakdown alpha (%g) < 0", (double)alpha);
180:         ksp->reason = KSP_DIVERGED_BREAKDOWN;
181:         break;
182:       }
183:     }
184:     if (alpha > 0.0) {
185:       if (!nopreconditioner) {
186:         alpha = PetscSqrtReal(alpha);
187:         PetscCall(VecScale(Z, 1.0 / alpha));
188:       }
189:       PetscCall(VecScale(V1, 1.0 / alpha)); /* alpha*V1 = Amat^T*U1 - beta*V */
190:     }
191:     rho    = PetscSqrtScalar(rhobar * rhobar + beta * beta);
192:     c      = rhobar / rho;
193:     s      = beta / rho;
194:     theta  = s * alpha;
195:     rhobar = -c * alpha;
196:     phi    = c * phibar;
197:     phibar = s * phibar;
198:     tau    = s * phi;

200:     PetscCall(VecAXPY(X, phi / rho, W)); /*    x <- x + (phi/rho) w   */

202:     if (lsqr->se) {
203:       PetscCall(VecCopy(W, W2));
204:       PetscCall(VecSquare(W2));
205:       PetscCall(VecScale(W2, 1.0 / (rho * rho)));
206:       PetscCall(VecAXPY(lsqr->se, 1.0, W2)); /* lsqr->se <- lsqr->se + (w^2/rho^2) */
207:     }
208:     if (nopreconditioner) {
209:       PetscCall(VecAYPX(W, -theta / rho, V1)); /* w <- v - (theta/rho) w */
210:     } else {
211:       PetscCall(VecAYPX(W, -theta / rho, Z)); /* w <- z - (theta/rho) w */
212:     }

214:     lsqr->arnorm = alpha * PetscAbsScalar(tau);
215:     rnorm        = PetscRealPart(phibar);

217:     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
218:     ksp->its++;
219:     ksp->rnorm = rnorm;
220:     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
221:     PetscCall(KSPLogResidualHistory(ksp, rnorm));
222:     PetscCall(KSPMonitor(ksp, i + 1, rnorm));
223:     PetscCall((*ksp->converged)(ksp, i + 1, rnorm, &ksp->reason, ksp->cnvP));
224:     if (ksp->reason) break;
225:     SWAP(U1, U, TMP);
226:     SWAP(V1, V, TMP);

228:     i++;
229:   } while (i < ksp->max_it);
230:   if (i >= ksp->max_it && !ksp->reason) ksp->reason = KSP_DIVERGED_ITS;

232:   /* Finish off the standard error estimates */
233:   if (lsqr->se) {
234:     tmp = 1.0;
235:     PetscCall(MatGetSize(Amat, &size1, &size2));
236:     if (size1 > size2) tmp = size1 - size2;
237:     tmp = rnorm / PetscSqrtScalar(tmp);
238:     PetscCall(VecSqrtAbs(lsqr->se));
239:     PetscCall(VecScale(lsqr->se, tmp));
240:   }
241:   PetscFunctionReturn(PETSC_SUCCESS);
242: }

244: static PetscErrorCode KSPDestroy_LSQR(KSP ksp)
245: {
246:   KSP_LSQR *lsqr = (KSP_LSQR *)ksp->data;

248:   PetscFunctionBegin;
249:   /* Free work vectors */
250:   if (lsqr->vwork_n) PetscCall(VecDestroyVecs(lsqr->nwork_n, &lsqr->vwork_n));
251:   if (lsqr->vwork_m) PetscCall(VecDestroyVecs(lsqr->nwork_m, &lsqr->vwork_m));
252:   PetscCall(VecDestroy(&lsqr->se));
253:   /* Revert convergence test */
254:   PetscCall(KSPSetConvergenceTest(ksp, lsqr->converged, lsqr->cnvP, lsqr->convergeddestroy));
255:   /* Free the KSP_LSQR context */
256:   PetscCall(PetscFree(ksp->data));
257:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPLSQRMonitorResidual_C", NULL));
258:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPLSQRMonitorResidualDrawLG_C", NULL));
259:   PetscFunctionReturn(PETSC_SUCCESS);
260: }

262: /*@
263:   KSPLSQRSetComputeStandardErrorVec - Compute a vector of standard error estimates during `KSPSolve()` for  `KSPLSQR`.

265:   Logically Collective

267:   Input Parameters:
268: + ksp - iterative context
269: - flg - compute the vector of standard estimates or not

271:   Level: intermediate

273:   Developer Notes:
274:   Vaclav: I'm not sure whether this vector is useful for anything.

276: .seealso: [](ch_ksp), `KSPSolve()`, `KSPLSQR`, `KSPLSQRGetStandardErrorVec()`
277: @*/
278: PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP ksp, PetscBool flg)
279: {
280:   KSP_LSQR *lsqr = (KSP_LSQR *)ksp->data;

282:   PetscFunctionBegin;
283:   lsqr->se_flg = flg;
284:   PetscFunctionReturn(PETSC_SUCCESS);
285: }

287: /*@
288:   KSPLSQRSetExactMatNorm - Compute exact matrix norm instead of iteratively refined estimate.

290:   Not Collective

292:   Input Parameters:
293: + ksp - iterative context
294: - flg - compute exact matrix norm or not

296:   Level: intermediate

298:   Notes:
299:   By default, `flg` = `PETSC_FALSE`. This is usually preferred to avoid possibly expensive computation of the norm.
300:   For `flg` = `PETSC_TRUE`, we call `MatNorm`(Amat,`NORM_FROBENIUS`,&lsqr->anorm) which will work only for some types of explicitly assembled matrices.
301:   This can affect convergence rate as `KSPLSQRConvergedDefault()` assumes different value of $||A||$ used in normal equation stopping criterion.

303: .seealso: [](ch_ksp), `KSPSolve()`, `KSPLSQR`, `KSPLSQRGetNorms()`, `KSPLSQRConvergedDefault()`
304: @*/
305: PetscErrorCode KSPLSQRSetExactMatNorm(KSP ksp, PetscBool flg)
306: {
307:   KSP_LSQR *lsqr = (KSP_LSQR *)ksp->data;

309:   PetscFunctionBegin;
310:   lsqr->exact_norm = flg;
311:   PetscFunctionReturn(PETSC_SUCCESS);
312: }

314: /*@
315:   KSPLSQRGetStandardErrorVec - Get vector of standard error estimates.
316:   Only available if -ksp_lsqr_set_standard_error was set to true
317:   or `KSPLSQRSetComputeStandardErrorVec`(ksp, `PETSC_TRUE`) was called.
318:   Otherwise returns `NULL`.

320:   Not Collective

322:   Input Parameter:
323: . ksp - iterative context

325:   Output Parameter:
326: . se - vector of standard estimates

328:   Level: intermediate

330:   Developer Notes:
331:   Vaclav: I'm not sure whether this vector is useful for anything.

333: .seealso: [](ch_ksp), `KSPSolve()`, `KSPLSQR`, `KSPLSQRSetComputeStandardErrorVec()`
334: @*/
335: PetscErrorCode KSPLSQRGetStandardErrorVec(KSP ksp, Vec *se)
336: {
337:   KSP_LSQR *lsqr = (KSP_LSQR *)ksp->data;

339:   PetscFunctionBegin;
340:   *se = lsqr->se;
341:   PetscFunctionReturn(PETSC_SUCCESS);
342: }

344: /*@
345:   KSPLSQRGetNorms - Get the norm estimates that `KSPLSQR` computes internally during `KSPSolve()`.

347:   Not Collective

349:   Input Parameter:
350: . ksp - iterative context

352:   Output Parameters:
353: + arnorm - good estimate of $\|(A*Pmat^{-T})*r\|$, where $r = A x - b$, used in specific stopping criterion
354: - anorm  - poor estimate of $\|A*Pmat^{-T}\|_{frobenius}$ used in specific stopping criterion

356:   Level: intermediate

358:   Notes:
359:   Output parameters are meaningful only after `KSPSolve()`.

361:   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.

363:   If `-ksp_lsqr_exact_mat_norm` is set or `KSPLSQRSetExactMatNorm`(ksp, `PETSC_TRUE`) called, then `anorm` is the exact Frobenius norm.

365: .seealso: [](ch_ksp), `KSPSolve()`, `KSPLSQR`, `KSPLSQRSetExactMatNorm()`
366: @*/
367: PetscErrorCode KSPLSQRGetNorms(KSP ksp, PetscReal *arnorm, PetscReal *anorm)
368: {
369:   KSP_LSQR *lsqr = (KSP_LSQR *)ksp->data;

371:   PetscFunctionBegin;
372:   if (arnorm) *arnorm = lsqr->arnorm;
373:   if (anorm) *anorm = lsqr->anorm;
374:   PetscFunctionReturn(PETSC_SUCCESS);
375: }

377: static PetscErrorCode KSPLSQRMonitorResidual_LSQR(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
378: {
379:   KSP_LSQR         *lsqr   = (KSP_LSQR *)ksp->data;
380:   PetscViewer       viewer = vf->viewer;
381:   PetscViewerFormat format = vf->format;
382:   char              normtype[256];
383:   PetscInt          tablevel;
384:   const char       *prefix;

386:   PetscFunctionBegin;
387:   PetscCall(PetscObjectGetTabLevel((PetscObject)ksp, &tablevel));
388:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ksp, &prefix));
389:   PetscCall(PetscStrncpy(normtype, KSPNormTypes[ksp->normtype], sizeof(normtype)));
390:   PetscCall(PetscStrtolower(normtype));
391:   PetscCall(PetscViewerPushFormat(viewer, format));
392:   PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
393:   if (n == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual norm, norm of normal equations, and matrix norm for %s solve.\n", prefix));
394:   if (!n) {
395:     PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " KSP resid norm %14.12e\n", n, (double)rnorm));
396:   } else {
397:     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));
398:   }
399:   PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
400:   PetscCall(PetscViewerPopFormat(viewer));
401:   PetscFunctionReturn(PETSC_SUCCESS);
402: }

404: /*@C
405:   KSPLSQRMonitorResidual - Prints the residual norm, as well as the normal equation residual norm, at each iteration of an iterative solver for the `KSPLSQR` solver

407:   Collective

409:   Input Parameters:
410: + ksp   - iterative context
411: . n     - iteration number
412: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
413: - vf    - The viewer context

415:   Options Database Key:
416: . -ksp_lsqr_monitor - Activates `KSPLSQRMonitorResidual()`

418:   Level: intermediate

420: .seealso: [](ch_ksp), `KSPLSQR`, `KSPMonitorSet()`, `KSPMonitorResidual()`, `KSPMonitorTrueResidualMaxNorm()`, `KSPLSQRMonitorResidualDrawLG()`
421: @*/
422: PetscErrorCode KSPLSQRMonitorResidual(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
423: {
424:   PetscFunctionBegin;
426:   PetscAssertPointer(vf, 4);
428:   PetscTryMethod(ksp, "KSPLSQRMonitorResidual_C", (KSP, PetscInt, PetscReal, PetscViewerAndFormat *), (ksp, n, rnorm, vf));
429:   PetscFunctionReturn(PETSC_SUCCESS);
430: }

432: static PetscErrorCode KSPLSQRMonitorResidualDrawLG_LSQR(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
433: {
434:   KSP_LSQR          *lsqr   = (KSP_LSQR *)ksp->data;
435:   PetscViewer        viewer = vf->viewer;
436:   PetscViewerFormat  format = vf->format;
437:   KSPConvergedReason reason;
438:   PetscReal          x[2], y[2];
439:   PetscDrawLG        lg;

441:   PetscFunctionBegin;
442:   PetscCall(PetscViewerPushFormat(viewer, format));
443:   PetscCall(PetscViewerDrawGetDrawLG(viewer, 0, &lg));
444:   if (!n) PetscCall(PetscDrawLGReset(lg));
445:   x[0] = (PetscReal)n;
446:   if (rnorm > 0.0) y[0] = PetscLog10Real(rnorm);
447:   else y[0] = -15.0;
448:   x[1] = (PetscReal)n;
449:   if (lsqr->arnorm > 0.0) y[1] = PetscLog10Real(lsqr->arnorm);
450:   else y[1] = -15.0;
451:   PetscCall(PetscDrawLGAddPoint(lg, x, y));
452:   PetscCall(KSPGetConvergedReason(ksp, &reason));
453:   if (n <= 20 || !(n % 5) || reason) {
454:     PetscCall(PetscDrawLGDraw(lg));
455:     PetscCall(PetscDrawLGSave(lg));
456:   }
457:   PetscCall(PetscViewerPopFormat(viewer));
458:   PetscFunctionReturn(PETSC_SUCCESS);
459: }

461: /*@C
462:   KSPLSQRMonitorResidualDrawLG - Plots the true residual norm at each iteration of an iterative solver for the `KSPLSQR` solver

464:   Collective

466:   Input Parameters:
467: + ksp   - iterative context
468: . n     - iteration number
469: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
470: - vf    - The viewer context

472:   Options Database Key:
473: . -ksp_lsqr_monitor draw::draw_lg - Activates `KSPMonitorTrueResidualDrawLG()`

475:   Level: intermediate

477: .seealso: [](ch_ksp), `KSPLSQR`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `KSPLSQRMonitorResidual()`, `KSPLSQRMonitorResidualDrawLGCreate()`
478: @*/
479: PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
480: {
481:   PetscFunctionBegin;
483:   PetscAssertPointer(vf, 4);
485:   PetscTryMethod(ksp, "KSPLSQRMonitorResidualDrawLG_C", (KSP, PetscInt, PetscReal, PetscViewerAndFormat *), (ksp, n, rnorm, vf));
486:   PetscFunctionReturn(PETSC_SUCCESS);
487: }

489: /*@C
490:   KSPLSQRMonitorResidualDrawLGCreate - Creates the line graph object for the `KSPLSQR` residual and normal equation residual norm

492:   Collective

494:   Input Parameters:
495: + viewer - The `PetscViewer`
496: . format - The viewer format
497: - ctx    - An optional user context

499:   Output Parameter:
500: . vf - The `PetscViewerAndFormat`

502:   Level: intermediate

504: .seealso: [](ch_ksp), `KSPLSQR`, `KSPMonitorSet()`, `KSPLSQRMonitorResidual()`, `KSPLSQRMonitorResidualDrawLG()`
505: @*/
506: PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer viewer, PetscViewerFormat format, void *ctx, PetscViewerAndFormat **vf)
507: {
508:   const char *names[] = {"residual", "normal eqn residual"};

510:   PetscFunctionBegin;
511:   PetscCall(PetscViewerAndFormatCreate(viewer, format, vf));
512:   (*vf)->data = ctx;
513:   PetscCall(PetscViewerMonitorLGSetUp(viewer, NULL, NULL, "Log Residual Norm", 2, names, PETSC_DECIDE, PETSC_DECIDE, 400, 300));
514:   PetscFunctionReturn(PETSC_SUCCESS);
515: }

517: static PetscErrorCode KSPSetFromOptions_LSQR(KSP ksp, PetscOptionItems PetscOptionsObject)
518: {
519:   KSP_LSQR *lsqr = (KSP_LSQR *)ksp->data;

521:   PetscFunctionBegin;
522:   PetscOptionsHeadBegin(PetscOptionsObject, "KSP LSQR Options");
523:   PetscCall(PetscOptionsBool("-ksp_lsqr_compute_standard_error", "Set Standard Error Estimates of Solution", "KSPLSQRSetComputeStandardErrorVec", lsqr->se_flg, &lsqr->se_flg, NULL));
524:   PetscCall(PetscOptionsBool("-ksp_lsqr_exact_mat_norm", "Compute exact matrix norm instead of iteratively refined estimate", "KSPLSQRSetExactMatNorm", lsqr->exact_norm, &lsqr->exact_norm, NULL));
525:   PetscCall(KSPMonitorSetFromOptions(ksp, "-ksp_lsqr_monitor", "lsqr_residual", NULL));
526:   PetscOptionsHeadEnd();
527:   PetscFunctionReturn(PETSC_SUCCESS);
528: }

530: static PetscErrorCode KSPView_LSQR(KSP ksp, PetscViewer viewer)
531: {
532:   KSP_LSQR *lsqr = (KSP_LSQR *)ksp->data;
533:   PetscBool isascii;

535:   PetscFunctionBegin;
536:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
537:   if (isascii) {
538:     if (lsqr->se) {
539:       PetscReal rnorm;
540:       PetscCall(VecNorm(lsqr->se, NORM_2, &rnorm));
541:       PetscCall(PetscViewerASCIIPrintf(viewer, "  norm of standard error %g, iterations %" PetscInt_FMT "\n", (double)rnorm, ksp->its));
542:     } else {
543:       PetscCall(PetscViewerASCIIPrintf(viewer, "  standard error not computed\n"));
544:     }
545:     if (lsqr->exact_norm) {
546:       PetscCall(PetscViewerASCIIPrintf(viewer, "  using exact matrix norm\n"));
547:     } else {
548:       PetscCall(PetscViewerASCIIPrintf(viewer, "  using inexact matrix norm\n"));
549:     }
550:   }
551:   PetscFunctionReturn(PETSC_SUCCESS);
552: }

554: /*@C
555:   KSPLSQRConvergedDefault - Determines convergence of the `KSPLSQR` Krylov method, including a check on the residual norm of the normal equations.

557:   Collective

559:   Input Parameters:
560: + ksp   - iterative context
561: . n     - iteration number
562: . rnorm - 2-norm residual value (may be estimated)
563: - ctx   - convergence context which must have been created by `KSPConvergedDefaultCreate()`

565:   Output Parameter:
566: . reason - the convergence reason

568:   Level: advanced

570:   Notes:
571:   This is not called directly but rather is passed to `KSPSetConvergenceTest()`. It is used automatically by `KSPLSQR`

573:   `KSPConvergedDefault()` is called first to check for convergence in $A*x=b$.
574:   If that does not determine convergence then checks convergence for the least squares problem, i.e., in $ \min_x |b - A x| $.
575:   Possible convergence for the least squares problem (which is based on the residual of the normal equations) are `KSP_CONVERGED_RTOL_NORMAL_EQUATIONS`
576:   and `KSP_CONVERGED_ATOL_NORMAL_EQUATIONS`.

578:   `KSP_CONVERGED_RTOL_NORMAL_EQUATIONS` is returned if $||A^T r|| < rtol ||A|| ||r||$.
579:   The matrix norm $||A||$ is an iteratively refined estimate, see `KSPLSQRGetNorms()`.
580:   This criterion is largely compatible with that in MATLAB `lsqr()`.

582: .seealso: [](ch_ksp), `KSPLSQR`, `KSPSetConvergenceTest()`, `KSPSetTolerances()`, `KSPConvergedSkip()`, `KSPConvergedReason`, `KSPGetConvergedReason()`,
583:           `KSPConvergedDefaultSetUIRNorm()`, `KSPConvergedDefaultSetUMIRNorm()`, `KSPConvergedDefaultCreate()`, `KSPConvergedDefaultDestroy()`,
584:           `KSPConvergedDefault()`, `KSPLSQRGetNorms()`, `KSPLSQRSetExactMatNorm()`
585: @*/
586: PetscErrorCode KSPLSQRConvergedDefault(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, void *ctx)
587: {
588:   KSP_LSQR *lsqr = (KSP_LSQR *)ksp->data;
589:   PetscReal xnorm;

591:   PetscFunctionBegin;
592:   /* check for convergence in A*x=b */
593:   PetscCall(KSPConvergedDefault(ksp, n, rnorm, reason, ctx));
594:   if (!n || *reason) PetscFunctionReturn(PETSC_SUCCESS);

596:   PetscCall(VecNorm(ksp->vec_sol, NORM_2, &xnorm));
597:   /* check for convergence in min{|b-A*x|} */
598:   if (lsqr->arnorm < ksp->rtol * ksp->rnorm0 + ksp->abstol * lsqr->anorm * xnorm) {
599:     PetscCall(PetscInfo(ksp, "LSQR solver has converged. Normal equation residual %14.12e is less than relative tolerance %14.12e times initial rhs norm %14.12e + absolute tolerance %14.12e times %s Frobenius norm of matrix %14.12e times solution %14.12e at iteration %" PetscInt_FMT "\n",
600:                         (double)lsqr->arnorm, (double)ksp->rtol, (double)ksp->rnorm0, (double)ksp->abstol, lsqr->exact_norm ? "exact" : "approx.", (double)lsqr->anorm, (double)xnorm, n));
601:     *reason = KSP_CONVERGED_RTOL_NORMAL_EQUATIONS;
602:   } else if (lsqr->arnorm < ksp->abstol * lsqr->anorm * rnorm) {
603:     PetscCall(PetscInfo(ksp, "LSQR solver has converged. Normal equation residual %14.12e is less than absolute tolerance %14.12e times %s Frobenius norm of matrix %14.12e times residual %14.12e at iteration %" PetscInt_FMT "\n", (double)lsqr->arnorm,
604:                         (double)ksp->abstol, lsqr->exact_norm ? "exact" : "approx.", (double)lsqr->anorm, (double)rnorm, n));
605:     *reason = KSP_CONVERGED_ATOL_NORMAL_EQUATIONS;
606:   }
607:   PetscFunctionReturn(PETSC_SUCCESS);
608: }

610: /*MC
611:    KSPLSQR - Implements LSQR  {cite}`paige.saunders:lsqr`

613:    Options Database Keys:
614: +  -ksp_lsqr_set_standard_error - set standard error estimates of solution, see `KSPLSQRSetComputeStandardErrorVec()` and `KSPLSQRGetStandardErrorVec()`
615: .  -ksp_lsqr_exact_mat_norm     - compute the exact matrix norm instead of using an iteratively refined estimate, see `KSPLSQRSetExactMatNorm()`
616: -  -ksp_lsqr_monitor            - monitor residual norm, norm of residual of normal equations $A^T A x = A^T b $, and estimate of matrix norm $||A||$

618:    Level: beginner

620:    Notes:
621:    Supports non-square (rectangular) matrices.  See `PETSCREGRESSORLINEAR` for the PETSc toolkit for solving linear regression problems, including least squares.

623:    This variant, when applied with no preconditioning is identical to the original published algorithm in exact arithmetic; however, in practice, with no preconditioning
624:    due to inexact arithmetic, it can converge differently. Hence when no preconditioner is used (`PCType` `PCNONE`) it automatically reverts to the original algorithm.

626:    With the PETSc built-in preconditioners, such as `PCICC`, one should call `KSPSetOperators`(ksp,A,A^T*A)) since the preconditioner needs to work
627:    for the normal equations $^T A$. For example, use `MatCreateNormal()`.

629:    Supports only left preconditioning.

631:    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()`.
632:    see `KSPLSQRConvergedDefault()`.

634:    In exact arithmetic the LSQR method (with no preconditioning) is identical to the `KSPCG` algorithm applied to the normal equations.
635:    The preconditioned variant was implemented by Bas van't Hof and is essentially a left preconditioning for the normal equations.
636:    It appears the implementation with preconditioning tracks the true (unpreconditioned) norm of the residual and uses that in the convergence test.

638:    Developer Note:
639:    How is this related to the `KSPCGNE` implementation? One difference is that `KSPCGNE` applies
640:    the preconditioner transpose times the preconditioner,  so one does not need to pass $A^T*A$ as the third argument to `KSPSetOperators()`.

642: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPSolve()`, `KSPLSQRConvergedDefault()`, `KSPLSQRSetComputeStandardErrorVec()`, `KSPLSQRGetStandardErrorVec()`, `KSPLSQRSetExactMatNorm()`, `KSPLSQRMonitorResidualDrawLGCreate()`, `KSPLSQRMonitorResidualDrawLG()`, `KSPLSQRMonitorResidual()`, `PETSCREGRESSORLINEAR`
643: M*/
644: PETSC_EXTERN PetscErrorCode KSPCreate_LSQR(KSP ksp)
645: {
646:   KSP_LSQR *lsqr;
647:   void     *ctx;

649:   PetscFunctionBegin;
650:   PetscCall(PetscNew(&lsqr));
651:   lsqr->se         = NULL;
652:   lsqr->se_flg     = PETSC_FALSE;
653:   lsqr->exact_norm = PETSC_FALSE;
654:   lsqr->anorm      = -1.0;
655:   lsqr->arnorm     = -1.0;
656:   ksp->data        = (void *)lsqr;
657:   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_LEFT, 3));

659:   ksp->ops->setup          = KSPSetUp_LSQR;
660:   ksp->ops->solve          = KSPSolve_LSQR;
661:   ksp->ops->destroy        = KSPDestroy_LSQR;
662:   ksp->ops->setfromoptions = KSPSetFromOptions_LSQR;
663:   ksp->ops->view           = KSPView_LSQR;

665:   /* Backup current convergence test; remove destroy routine from KSP to prevent destroying the convergence context in KSPSetConvergenceTest() */
666:   PetscCall(KSPGetAndClearConvergenceTest(ksp, &lsqr->converged, &lsqr->cnvP, &lsqr->convergeddestroy));
667:   /* Override current convergence test */
668:   PetscCall(KSPConvergedDefaultCreate(&ctx));
669:   PetscCall(KSPSetConvergenceTest(ksp, KSPLSQRConvergedDefault, ctx, KSPConvergedDefaultDestroy));
670:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPLSQRMonitorResidual_C", KSPLSQRMonitorResidual_LSQR));
671:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPLSQRMonitorResidualDrawLG_C", KSPLSQRMonitorResidualDrawLG_LSQR));
672:   PetscFunctionReturn(PETSC_SUCCESS);
673: }