Actual source code: lsqr.c


  2: /* lourens.vanzanen@shell.com contributed the standard error estimates of the solution, Jul 25, 2006 */
  3: /* Bas van't Hof contributed the preconditioned aspects Feb 10, 2010 */

  5: #define SWAP(a,b,c) { c = a; a = b; b = c; }

  7: #include <petsc/private/kspimpl.h>
  8: #include <petscdraw.h>

 10: typedef struct {
 11:   PetscInt  nwork_n,nwork_m;
 12:   Vec       *vwork_m;   /* work vectors of length m, where the system is size m x n */
 13:   Vec       *vwork_n;   /* work vectors of length n */
 14:   Vec       se;         /* Optional standard error vector */
 15:   PetscBool se_flg;     /* flag for -ksp_lsqr_set_standard_error */
 16:   PetscBool exact_norm; /* flag for -ksp_lsqr_exact_mat_norm */
 17:   PetscReal arnorm;     /* Good estimate of norm((A*inv(Pmat))'*r), where r = A*x - b, used in specific stopping criterion */
 18:   PetscReal anorm;      /* Poor estimate of norm(A*inv(Pmat),'fro') used in specific stopping criterion */
 19:   /* Backup previous convergence test */
 20:   PetscErrorCode        (*converged)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
 21:   PetscErrorCode        (*convergeddestroy)(void*);
 22:   void                  *cnvP;
 23: } KSP_LSQR;

 25: static PetscErrorCode  VecSquare(Vec v)
 26: {
 28:   PetscScalar    *x;
 29:   PetscInt       i, n;

 32:   VecGetLocalSize(v, &n);
 33:   VecGetArray(v, &x);
 34:   for (i = 0; i < n; i++) x[i] *= PetscConj(x[i]);
 35:   VecRestoreArray(v, &x);
 36:   return(0);
 37: }

 39: static PetscErrorCode KSPSetUp_LSQR(KSP ksp)
 40: {
 42:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;
 43:   PetscBool      nopreconditioner;

 46:   PetscObjectTypeCompare((PetscObject)ksp->pc,PCNONE,&nopreconditioner);

 48:   if (lsqr->vwork_m) {
 49:     VecDestroyVecs(lsqr->nwork_m,&lsqr->vwork_m);
 50:   }

 52:   if (lsqr->vwork_n) {
 53:     VecDestroyVecs(lsqr->nwork_n,&lsqr->vwork_n);
 54:   }

 56:   lsqr->nwork_m = 2;
 57:   if (nopreconditioner) lsqr->nwork_n = 4;
 58:   else lsqr->nwork_n = 5;
 59:   KSPCreateVecs(ksp,lsqr->nwork_n,&lsqr->vwork_n,lsqr->nwork_m,&lsqr->vwork_m);

 61:   if (lsqr->se_flg && !lsqr->se) {
 62:     VecDuplicate(lsqr->vwork_n[0],&lsqr->se);
 63:     VecSet(lsqr->se,PETSC_INFINITY);
 64:   } else if (!lsqr->se_flg) {
 65:     VecDestroy(&lsqr->se);
 66:   }
 67:   return(0);
 68: }

 70: static PetscErrorCode KSPSolve_LSQR(KSP ksp)
 71: {
 73:   PetscInt       i,size1,size2;
 74:   PetscScalar    rho,rhobar,phi,phibar,theta,c,s,tmp,tau;
 75:   PetscReal      beta,alpha,rnorm;
 76:   Vec            X,B,V,V1,U,U1,TMP,W,W2,Z = NULL;
 77:   Mat            Amat,Pmat;
 78:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;
 79:   PetscBool      diagonalscale,nopreconditioner;

 82:   PCGetDiagonalScale(ksp->pc,&diagonalscale);
 83:   if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);

 85:   PCGetOperators(ksp->pc,&Amat,&Pmat);
 86:   PetscObjectTypeCompare((PetscObject)ksp->pc,PCNONE,&nopreconditioner);

 88:   /* vectors of length m, where system size is mxn */
 89:   B  = ksp->vec_rhs;
 90:   U  = lsqr->vwork_m[0];
 91:   U1 = lsqr->vwork_m[1];

 93:   /* vectors of length n */
 94:   X  = ksp->vec_sol;
 95:   W  = lsqr->vwork_n[0];
 96:   V  = lsqr->vwork_n[1];
 97:   V1 = lsqr->vwork_n[2];
 98:   W2 = lsqr->vwork_n[3];
 99:   if (!nopreconditioner) Z = lsqr->vwork_n[4];

101:   /* standard error vector */
102:   if (lsqr->se) {
103:     VecSet(lsqr->se,0.0);
104:   }

106:   /* Compute initial residual, temporarily use work vector u */
107:   if (!ksp->guess_zero) {
108:     KSP_MatMult(ksp,Amat,X,U);       /*   u <- b - Ax     */
109:     VecAYPX(U,-1.0,B);
110:   } else {
111:     VecCopy(B,U);            /*   u <- b (x is 0) */
112:   }

114:   /* Test for nothing to do */
115:   VecNorm(U,NORM_2,&rnorm);
116:   KSPCheckNorm(ksp,rnorm);
117:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
118:   ksp->its   = 0;
119:   ksp->rnorm = rnorm;
120:   PetscObjectSAWsGrantAccess((PetscObject)ksp);
121:   KSPLogResidualHistory(ksp,rnorm);
122:   KSPMonitor(ksp,0,rnorm);
123:   (*ksp->converged)(ksp,0,rnorm,&ksp->reason,ksp->cnvP);
124:   if (ksp->reason) return(0);

126:   beta = rnorm;
127:   VecScale(U,1.0/beta);
128:   KSP_MatMultTranspose(ksp,Amat,U,V);
129:   if (nopreconditioner) {
130:     VecNorm(V,NORM_2,&alpha);
131:     KSPCheckNorm(ksp,rnorm);
132:   } else {
133:     /* this is an application of the preconditioner for the normal equations; not the operator, see the manual page */
134:     PCApply(ksp->pc,V,Z);
135:     VecDotRealPart(V,Z,&alpha);
136:     if (alpha <= 0.0) {
137:       ksp->reason = KSP_DIVERGED_BREAKDOWN;
138:       return(0);
139:     }
140:     alpha = PetscSqrtReal(alpha);
141:     VecScale(Z,1.0/alpha);
142:   }
143:   VecScale(V,1.0/alpha);

145:   if (nopreconditioner) {
146:     VecCopy(V,W);
147:   } else {
148:     VecCopy(Z,W);
149:   }

151:   if (lsqr->exact_norm) {
152:     MatNorm(Amat,NORM_FROBENIUS,&lsqr->anorm);
153:   } else lsqr->anorm = 0.0;

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

175:     KSP_MatMultTranspose(ksp,Amat,U1,V1);
176:     VecAXPY(V1,-beta,V);
177:     if (nopreconditioner) {
178:       VecNorm(V1,NORM_2,&alpha);
179:       KSPCheckNorm(ksp,alpha);
180:     } else {
181:       PCApply(ksp->pc,V1,Z);
182:       VecDotRealPart(V1,Z,&alpha);
183:       if (alpha <= 0.0) {
184:         ksp->reason = KSP_DIVERGED_BREAKDOWN;
185:         break;
186:       }
187:       alpha = PetscSqrtReal(alpha);
188:       VecScale(Z,1.0/alpha);
189:     }
190:     VecScale(V1,1.0/alpha); /* alpha*V1 = Amat^T*U1 - beta*V */
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:     VecAXPY(X,phi/rho,W);  /*    x <- x + (phi/rho) w   */

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

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

217:     PetscObjectSAWsTakeAccess((PetscObject)ksp);
218:     ksp->its++;
219:     ksp->rnorm = rnorm;
220:     PetscObjectSAWsGrantAccess((PetscObject)ksp);
221:     KSPLogResidualHistory(ksp,rnorm);
222:     KSPMonitor(ksp,i+1,rnorm);
223:     (*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:     MatGetSize(Amat,&size1,&size2);
236:     if (size1 > size2) tmp = size1 - size2;
237:     tmp  = rnorm / PetscSqrtScalar(tmp);
238:     VecSqrtAbs(lsqr->se);
239:     VecScale(lsqr->se,tmp);
240:   }
241:   return(0);
242: }


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

251:   /* Free work vectors */
252:   if (lsqr->vwork_n) {
253:     VecDestroyVecs(lsqr->nwork_n,&lsqr->vwork_n);
254:   }
255:   if (lsqr->vwork_m) {
256:     VecDestroyVecs(lsqr->nwork_m,&lsqr->vwork_m);
257:   }
258:   VecDestroy(&lsqr->se);
259:   /* Revert convergence test */
260:   KSPSetConvergenceTest(ksp,lsqr->converged,lsqr->cnvP,lsqr->convergeddestroy);
261:   /* Free the KSP_LSQR context */
262:   PetscFree(ksp->data);
263:   PetscObjectComposeFunction((PetscObject)ksp,"KSPLSQRMonitorResidual_C",NULL);
264:   PetscObjectComposeFunction((PetscObject)ksp,"KSPLSQRMonitorResidualDrawLG_C",NULL);
265:   return(0);
266: }

268: /*@
269:    KSPLSQRSetComputeStandardErrorVec - Compute vector of standard error estimates during KSPSolve_LSQR().

271:    Not Collective

273:    Input Parameters:
274: +  ksp   - iterative context
275: -  flg   - compute the vector of standard estimates or not

277:    Developer notes:
278:    Vaclav: I'm not sure whether this vector is useful for anything.

280:    Level: intermediate

282: .seealso: KSPSolve(), KSPLSQR, KSPLSQRGetStandardErrorVec()
283: @*/
284: PetscErrorCode  KSPLSQRSetComputeStandardErrorVec(KSP ksp, PetscBool flg)
285: {
286:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;

289:   lsqr->se_flg = flg;
290:   return(0);
291: }

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

296:    Not Collective

298:    Input Parameters:
299: +  ksp   - iterative context
300: -  flg   - compute exact matrix norm or not

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

307:    Level: intermediate

309: .seealso: KSPSolve(), KSPLSQR, KSPLSQRGetNorms(), KSPLSQRConvergedDefault()
310: @*/
311: PetscErrorCode  KSPLSQRSetExactMatNorm(KSP ksp, PetscBool flg)
312: {
313:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;

316:   lsqr->exact_norm = flg;
317:   return(0);
318: }

320: /*@
321:    KSPLSQRGetStandardErrorVec - Get vector of standard error estimates.
322:    Only available if -ksp_lsqr_set_standard_error was set to true
323:    or KSPLSQRSetComputeStandardErrorVec(ksp, PETSC_TRUE) was called.
324:    Otherwise returns NULL.

326:    Not Collective

328:    Input Parameters:
329: .  ksp   - iterative context

331:    Output Parameters:
332: .  se - vector of standard estimates

334:    Options Database Keys:
335: .   -ksp_lsqr_set_standard_error  - set standard error estimates of solution

337:    Developer notes:
338:    Vaclav: I'm not sure whether this vector is useful for anything.

340:    Level: intermediate

342: .seealso: KSPSolve(), KSPLSQR, KSPLSQRSetComputeStandardErrorVec()
343: @*/
344: PetscErrorCode  KSPLSQRGetStandardErrorVec(KSP ksp,Vec *se)
345: {
346:   KSP_LSQR *lsqr = (KSP_LSQR*)ksp->data;

349:   *se = lsqr->se;
350:   return(0);
351: }

353: /*@
354:    KSPLSQRGetNorms - Get norm estimates that LSQR computes internally during KSPSolve().

356:    Not Collective

358:    Input Parameters:
359: .  ksp   - iterative context

361:    Output Parameters:
362: +  arnorm - good estimate of norm((A*inv(Pmat))'*r), where r = A*x - b, used in specific stopping criterion
363: -  anorm - poor estimate of norm(A*inv(Pmat),'fro') used in specific stopping criterion

365:    Notes:
366:    Output parameters are meaningful only after KSPSolve().
367:    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.
368:    If -ksp_lsqr_exact_mat_norm is set or KSPLSQRSetExactMatNorm(ksp, PETSC_TRUE) called, then anorm is exact Frobenius norm.

370:    Level: intermediate

372: .seealso: KSPSolve(), KSPLSQR, KSPLSQRSetExactMatNorm()
373: @*/
374: PetscErrorCode  KSPLSQRGetNorms(KSP ksp,PetscReal *arnorm, PetscReal *anorm)
375: {
376:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;

379:   if (arnorm)   *arnorm = lsqr->arnorm;
380:   if (anorm)    *anorm = lsqr->anorm;
381:   return(0);
382: }

384: PetscErrorCode KSPLSQRMonitorResidual_LSQR(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
385: {
386:   KSP_LSQR          *lsqr  = (KSP_LSQR*)ksp->data;
387:   PetscViewer       viewer = vf->viewer;
388:   PetscViewerFormat format = vf->format;
389:   char              normtype[256];
390:   PetscInt          tablevel;
391:   const char        *prefix;
392:   PetscErrorCode    ierr;

395:   PetscObjectGetTabLevel((PetscObject) ksp, &tablevel);
396:   PetscObjectGetOptionsPrefix((PetscObject) ksp, &prefix);
397:   PetscStrncpy(normtype, KSPNormTypes[ksp->normtype], sizeof(normtype));
398:   PetscStrtolower(normtype);
399:   PetscViewerPushFormat(viewer, format);
400:   PetscViewerASCIIAddTab(viewer, tablevel);
401:   if (n == 0 && prefix) {PetscViewerASCIIPrintf(viewer, "  Residual norm, norm of normal equations, and matrix norm for %s solve.\n", prefix);}
402:   if (!n) {
403:     PetscViewerASCIIPrintf(viewer,"%3D KSP resid norm %14.12e\n",n,(double)rnorm);
404:   } else {
405:     PetscViewerASCIIPrintf(viewer,"%3D 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);
406:   }
407:   PetscViewerASCIISubtractTab(viewer, tablevel);
408:   PetscViewerPopFormat(viewer);
409:   return(0);
410: }

412: /*@C
413:   KSPLSQRMonitorResidual - Prints the residual norm, as well as the normal equation residual norm, at each iteration of an iterative solver.

415:   Collective on ksp

417:   Input Parameters:
418: + ksp   - iterative context
419: . n     - iteration number
420: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
421: - vf    - The viewer context

423:   Options Database Key:
424: . -ksp_lsqr_monitor - Activates KSPLSQRMonitorResidual()

426:   Level: intermediate

428: .seealso: KSPMonitorSet(), KSPMonitorResidual(),KSPMonitorTrueResidualMaxNorm()
429: @*/
430: PetscErrorCode KSPLSQRMonitorResidual(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
431: {

438:   PetscTryMethod(ksp, "KSPLSQRMonitorResidual_C", (KSP,PetscInt,PetscReal,PetscViewerAndFormat*), (ksp,n,rnorm,vf));
439:   return(0);
440: }

442: PetscErrorCode KSPLSQRMonitorResidualDrawLG_LSQR(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
443: {
444:   KSP_LSQR           *lsqr  = (KSP_LSQR*)ksp->data;
445:   PetscViewer        viewer = vf->viewer;
446:   PetscViewerFormat  format = vf->format;
447:   PetscDrawLG        lg     = vf->lg;
448:   KSPConvergedReason reason;
449:   PetscReal          x[2], y[2];
450:   PetscErrorCode     ierr;

453:   PetscViewerPushFormat(viewer, format);
454:   if (!n) {PetscDrawLGReset(lg);}
455:   x[0] = (PetscReal) n;
456:   if (rnorm > 0.0) y[0] = PetscLog10Real(rnorm);
457:   else y[0] = -15.0;
458:   x[1] = (PetscReal) n;
459:   if (lsqr->arnorm > 0.0) y[1] = PetscLog10Real(lsqr->arnorm);
460:   else y[1] = -15.0;
461:   PetscDrawLGAddPoint(lg, x, y);
462:   KSPGetConvergedReason(ksp, &reason);
463:   if (n <= 20 || !(n % 5) || reason) {
464:     PetscDrawLGDraw(lg);
465:     PetscDrawLGSave(lg);
466:   }
467:   PetscViewerPopFormat(viewer);
468:   return(0);
469: }

471: /*@C
472:   KSPLSQRMonitorResidualDrawLG - Plots the true residual norm at each iteration of an iterative solver.

474:   Collective on ksp

476:   Input Parameters:
477: + ksp   - iterative context
478: . n     - iteration number
479: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
480: - vf    - The viewer context

482:   Options Database Key:
483: . -ksp_lsqr_monitor draw::draw_lg - Activates KSPMonitorTrueResidualDrawLG()

485:   Level: intermediate

487: .seealso: KSPMonitorSet(), KSPMonitorTrueResidual()
488: @*/
489: PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
490: {

498:   PetscTryMethod(ksp, "KSPLSQRMonitorResidualDrawLG_C", (KSP,PetscInt,PetscReal,PetscViewerAndFormat*), (ksp,n,rnorm,vf));
499:   return(0);
500: }

502: /*@C
503:   KSPLSQRMonitorResidualDrawLGCreate - Creates the plotter for the LSQR residual and normal eqn residual.

505:   Collective on ksp

507:   Input Parameters:
508: + viewer - The PetscViewer
509: . format - The viewer format
510: - ctx    - An optional user context

512:   Output Parameter:
513: . vf    - The viewer context

515:   Level: intermediate

517: .seealso: KSPMonitorSet(), KSPLSQRMonitorResidual()
518: @*/
519: PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer viewer, PetscViewerFormat format, void *ctx, PetscViewerAndFormat **vf)
520: {
521:   const char    *names[] = {"residual", "normal eqn residual"};

525:   PetscViewerAndFormatCreate(viewer, format, vf);
526:   (*vf)->data = ctx;
527:   KSPMonitorLGCreate(PetscObjectComm((PetscObject) viewer), NULL, NULL, "Log Residual Norm", 2, names, PETSC_DECIDE, PETSC_DECIDE, 400, 300, &(*vf)->lg);
528:   return(0);
529: }

531: PetscErrorCode KSPSetFromOptions_LSQR(PetscOptionItems *PetscOptionsObject,KSP ksp)
532: {
534:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;

537:   PetscOptionsHead(PetscOptionsObject,"KSP LSQR Options");
538:   PetscOptionsBool("-ksp_lsqr_compute_standard_error","Set Standard Error Estimates of Solution","KSPLSQRSetComputeStandardErrorVec",lsqr->se_flg,&lsqr->se_flg,NULL);
539:   PetscOptionsBool("-ksp_lsqr_exact_mat_norm","Compute exact matrix norm instead of iteratively refined estimate","KSPLSQRSetExactMatNorm",lsqr->exact_norm,&lsqr->exact_norm,NULL);
540:   KSPMonitorSetFromOptions(ksp, "-ksp_lsqr_monitor", "lsqr_residual", NULL);
541:   PetscOptionsTail();
542:   return(0);
543: }

545: PetscErrorCode KSPView_LSQR(KSP ksp,PetscViewer viewer)
546: {
547:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;
549:   PetscBool      iascii;

552:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
553:   if (iascii) {
554:     if (lsqr->se) {
555:       PetscReal rnorm;
556:       VecNorm(lsqr->se,NORM_2,&rnorm);
557:       PetscViewerASCIIPrintf(viewer,"  norm of standard error %g, iterations %d\n",(double)rnorm,ksp->its);
558:     } else {
559:       PetscViewerASCIIPrintf(viewer,"  standard error not computed\n");
560:     }
561:     if (lsqr->exact_norm) {
562:       PetscViewerASCIIPrintf(viewer,"  using exact matrix norm\n");
563:     } else {
564:       PetscViewerASCIIPrintf(viewer,"  using inexact matrix norm\n");
565:     }
566:   }
567:   return(0);
568: }

570: /*@C
571:    KSPLSQRConvergedDefault - Determines convergence of the LSQR Krylov method.

573:    Collective on ksp

575:    Input Parameters:
576: +  ksp   - iterative context
577: .  n     - iteration number
578: .  rnorm - 2-norm residual value (may be estimated)
579: -  ctx - convergence context which must be created by KSPConvergedDefaultCreate()

581:    reason is set to:
582: +   positive - if the iteration has converged;
583: .   negative - if residual norm exceeds divergence threshold;
584: -   0 - otherwise.

586:    Notes:
587:    KSPConvergedDefault() is called first to check for convergence in A*x=b.
588:    If that does not determine convergence then checks convergence for the least squares problem, i.e. in min{|b-A*x|}.
589:    Possible convergence for the least squares problem (which is based on the residual of the normal equations) are KSP_CONVERGED_RTOL_NORMAL norm and KSP_CONVERGED_ATOL_NORMAL.
590:    KSP_CONVERGED_RTOL_NORMAL is returned if ||A'*r|| < rtol * ||A|| * ||r||.
591:    Matrix norm ||A|| is iteratively refined estimate, see KSPLSQRGetNorms().
592:    This criterion is now largely compatible with that in MATLAB lsqr().

594:    Level: intermediate

596: .seealso: KSPLSQR, KSPSetConvergenceTest(), KSPSetTolerances(), KSPConvergedSkip(), KSPConvergedReason, KSPGetConvergedReason(),
597:           KSPConvergedDefaultSetUIRNorm(), KSPConvergedDefaultSetUMIRNorm(), KSPConvergedDefaultCreate(), KSPConvergedDefaultDestroy(), KSPConvergedDefault(), KSPLSQRGetNorms(), KSPLSQRSetExactMatNorm()
598: @*/
599: PetscErrorCode  KSPLSQRConvergedDefault(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *ctx)
600: {
602:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;

605:   /* check for convergence in A*x=b */
606:   KSPConvergedDefault(ksp,n,rnorm,reason,ctx);
607:   if (!n || *reason) return(0);

609:   /* check for convergence in min{|b-A*x|} */
610:   if (lsqr->arnorm < ksp->abstol) {
611:     PetscInfo3(ksp,"LSQR solver has converged. Normal equation residual %14.12e is less than absolute tolerance %14.12e at iteration %D\n",(double)lsqr->arnorm,(double)ksp->abstol,n);
612:     *reason = KSP_CONVERGED_ATOL_NORMAL;
613:   } else if (lsqr->arnorm < ksp->rtol * lsqr->anorm * rnorm) {
614:     PetscInfo6(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 %D\n",(double)lsqr->arnorm,(double)ksp->rtol,lsqr->exact_norm?"exact":"approx.",(double)lsqr->anorm,(double)rnorm,n);
615:     *reason = KSP_CONVERGED_RTOL_NORMAL;
616:   }
617:   return(0);
618: }

620: /*MC
621:      KSPLSQR - This implements LSQR

623:    Options Database Keys:
624: +   -ksp_lsqr_set_standard_error  - set standard error estimates of solution, see KSPLSQRSetComputeStandardErrorVec() and KSPLSQRGetStandardErrorVec()
625: .   -ksp_lsqr_exact_mat_norm - compute exact matrix norm instead of iteratively refined estimate, see KSPLSQRSetExactMatNorm()
626: -   -ksp_lsqr_monitor - monitor residual norm, norm of residual of normal equations A'*A x = A' b, and estimate of matrix norm ||A||

628:    Level: beginner

630:    Notes:
631:      Supports non-square (rectangular) matrices.

633:      This varient, when applied with no preconditioning is identical to the original algorithm in exact arithematic; however, in practice, with no preconditioning
634:      due to inexact arithematic, it can converge differently. Hence when no preconditioner is used (PCType PCNONE) it automatically reverts to the original algorithm.

636:      With the PETSc built-in preconditioners, such as ICC, one should call KSPSetOperators(ksp,A,A'*A)) since the preconditioner needs to work
637:      for the normal equations A'*A.

639:      Supports only left preconditioning.

641:      For least squares problems wit nonzero residual A*x - b, there are additional convergence tests for the residual of the normal equations, A'*(b - Ax), see KSPLSQRConvergedDefault().

643:    References:
644: .  1. - The original unpreconditioned algorithm can be found in Paige and Saunders, ACM Transactions on Mathematical Software, Vol 8, 1982.

646:      In exact arithmetic the LSQR method (with no preconditioning) is identical to the KSPCG algorithm applied to the normal equations.
647:      The preconditioned variant was implemented by Bas van't Hof and is essentially a left preconditioning for the Normal Equations. It appears the implementation with preconditioner
648:      track the true norm of the residual and uses that in the convergence test.

650:    Developer Notes:
651:     How is this related to the KSPCGNE implementation? One difference is that KSPCGNE applies
652:             the preconditioner transpose times the preconditioner,  so one does not need to pass A'*A as the third argument to KSPSetOperators().



656: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPSolve(), KSPLSQRConvergedDefault(), KSPLSQRSetComputeStandardErrorVec(), KSPLSQRGetStandardErrorVec(), KSPLSQRSetExactMatNorm()

658: M*/
659: PETSC_EXTERN PetscErrorCode KSPCreate_LSQR(KSP ksp)
660: {
661:   KSP_LSQR       *lsqr;
662:   void           *ctx;

666:   PetscNewLog(ksp,&lsqr);
667:   lsqr->se     = NULL;
668:   lsqr->se_flg = PETSC_FALSE;
669:   lsqr->exact_norm = PETSC_FALSE;
670:   lsqr->anorm  = -1.0;
671:   lsqr->arnorm = -1.0;
672:   ksp->data    = (void*)lsqr;
673:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,3);

675:   ksp->ops->setup          = KSPSetUp_LSQR;
676:   ksp->ops->solve          = KSPSolve_LSQR;
677:   ksp->ops->destroy        = KSPDestroy_LSQR;
678:   ksp->ops->setfromoptions = KSPSetFromOptions_LSQR;
679:   ksp->ops->view           = KSPView_LSQR;

681:   /* Backup current convergence test; remove destroy routine from KSP to prevent destroying the convergence context in KSPSetConvergenceTest() */
682:   KSPGetAndClearConvergenceTest(ksp,&lsqr->converged,&lsqr->cnvP,&lsqr->convergeddestroy);
683:   /* Override current convergence test */
684:   KSPConvergedDefaultCreate(&ctx);
685:   KSPSetConvergenceTest(ksp,KSPLSQRConvergedDefault,ctx,KSPConvergedDefaultDestroy);
686:   PetscObjectComposeFunction((PetscObject)ksp,"KSPLSQRMonitorResidual_C",KSPLSQRMonitorResidual_LSQR);
687:   PetscObjectComposeFunction((PetscObject)ksp,"KSPLSQRMonitorResidualDrawLG_C",KSPLSQRMonitorResidualDrawLG_LSQR);
688:   return(0);
689: }