Actual source code: lsqr.c

petsc-master 2016-09-25
Report Typos and Errors
  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>

  9: typedef struct {
 10:   PetscInt  nwork_n,nwork_m;
 11:   Vec       *vwork_m;   /* work vectors of length m, where the system is size m x n */
 12:   Vec       *vwork_n;   /* work vectors of length n */
 13:   Vec       se;         /* Optional standard error vector */
 14:   PetscBool se_flg;     /* flag for -ksp_lsqr_set_standard_error */
 15:   PetscReal arnorm;     /* Norm of the vector A.r */
 16:   PetscReal anorm;      /* Frobenius norm of the matrix A */
 17:   PetscReal rhs_norm;   /* Norm of the right hand side */
 18: } KSP_LSQR;

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

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

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

 45:   PetscObjectTypeCompare((PetscObject)ksp->pc,PCNONE,&nopreconditioner);
 46:   /*  nopreconditioner =PETSC_FALSE; */

 48:   lsqr->nwork_m = 2;
 49:   if (lsqr->vwork_m) {
 50:     VecDestroyVecs(lsqr->nwork_m,&lsqr->vwork_m);
 51:   }
 52:   if (nopreconditioner) lsqr->nwork_n = 4;
 53:   else lsqr->nwork_n = 5;

 55:   if (lsqr->vwork_n) {
 56:     VecDestroyVecs(lsqr->nwork_n,&lsqr->vwork_n);
 57:   }
 58:   KSPCreateVecs(ksp,lsqr->nwork_n,&lsqr->vwork_n,lsqr->nwork_m,&lsqr->vwork_m);
 59:   if (lsqr->se_flg && !lsqr->se) {
 60:     /* lsqr->se is not set by user, get it from pmat */
 61:     Vec *se;
 62:     KSPCreateVecs(ksp,1,&se,0,NULL);
 63:     lsqr->se = *se;
 64:     PetscFree(se);
 65:   }
 66:   return(0);
 67: }

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

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

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

 89:   /*  nopreconditioner =PETSC_FALSE; */
 90:   /* Calculate norm of right hand side */
 91:   VecNorm(ksp->vec_rhs,NORM_2,&lsqr->rhs_norm);

 93:   /* mark norm of matrix with negative number to indicate it has not yet been computed */
 94:   lsqr->anorm = -1.0;

 96:   /* vectors of length m, where system size is mxn */
 97:   B  = ksp->vec_rhs;
 98:   U  = lsqr->vwork_m[0];
 99:   U1 = lsqr->vwork_m[1];

101:   /* vectors of length n */
102:   X  = ksp->vec_sol;
103:   W  = lsqr->vwork_n[0];
104:   V  = lsqr->vwork_n[1];
105:   V1 = lsqr->vwork_n[2];
106:   W2 = lsqr->vwork_n[3];
107:   if (!nopreconditioner) Z = lsqr->vwork_n[4];

109:   /* standard error vector */
110:   SE = lsqr->se;
111:   if (SE) {
112:     VecGetSize(SE,&size1);
113:     VecGetSize(X,&size2);
114:     if (size1 != size2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Standard error vector (size %d) does not match solution vector (size %d)",size1,size2);
115:     VecSet(SE,0.0);
116:   }

118:   /* Compute initial residual, temporarily use work vector u */
119:   if (!ksp->guess_zero) {
120:     KSP_MatMult(ksp,Amat,X,U);       /*   u <- b - Ax     */
121:     VecAYPX(U,-1.0,B);
122:   } else {
123:     VecCopy(B,U);            /*   u <- b (x is 0) */
124:   }

126:   /* Test for nothing to do */
127:   VecNorm(U,NORM_2,&rnorm);
128:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
129:   ksp->its   = 0;
130:   ksp->rnorm = rnorm;
131:   PetscObjectSAWsGrantAccess((PetscObject)ksp);
132:   KSPLogResidualHistory(ksp,rnorm);
133:   KSPMonitor(ksp,0,rnorm);
134:   (*ksp->converged)(ksp,0,rnorm,&ksp->reason,ksp->cnvP);
135:   if (ksp->reason) return(0);

137:   beta = rnorm;
138:   VecScale(U,1.0/beta);
139:   KSP_MatMultTranspose(ksp,Amat,U,V);
140:   if (nopreconditioner) {
141:     VecNorm(V,NORM_2,&alpha);
142:   } else {
143:     PCApply(ksp->pc,V,Z);
144:     VecDotRealPart(V,Z,&alpha);
145:     if (alpha <= 0.0) {
146:       ksp->reason = KSP_DIVERGED_BREAKDOWN;
147:       return(0);
148:     }
149:     alpha = PetscSqrtReal(alpha);
150:     VecScale(Z,1.0/alpha);
151:   }
152:   VecScale(V,1.0/alpha);

154:   if (nopreconditioner) {
155:     VecCopy(V,W);
156:   } else {
157:     VecCopy(Z,W);
158:   }

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

176:     KSP_MatMultTranspose(ksp,Amat,U1,V1);
177:     VecAXPY(V1,-beta,V);
178:     if (nopreconditioner) {
179:       VecNorm(V1,NORM_2,&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 (SE) {
203:       VecCopy(W,W2);
204:       VecSquare(W2);
205:       VecScale(W2,1.0/(rho*rho));
206:       VecAXPY(SE, 1.0, W2); /* SE <- 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 (SE) {
234:     tmp  = 1.0;
235:     MatGetSize(Amat,&size1,&size2);
236:     if (size1 > size2) tmp = size1 - size2;
237:     tmp  = rnorm / PetscSqrtScalar(tmp);
238:     VecSqrtAbs(SE);
239:     VecScale(SE,tmp);
240:   }
241:   return(0);
242: }


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

253:   /* Free work vectors */
254:   if (lsqr->vwork_n) {
255:     VecDestroyVecs(lsqr->nwork_n,&lsqr->vwork_n);
256:   }
257:   if (lsqr->vwork_m) {
258:     VecDestroyVecs(lsqr->nwork_m,&lsqr->vwork_m);
259:   }
260:   if (lsqr->se_flg) {
261:     VecDestroy(&lsqr->se);
262:   }
263:   PetscFree(ksp->data);
264:   return(0);
265: }

269: PetscErrorCode  KSPLSQRSetStandardErrorVec(KSP ksp, Vec se)
270: {
271:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;

275:   VecDestroy(&lsqr->se);
276:   lsqr->se = se;
277:   return(0);
278: }

282: PetscErrorCode  KSPLSQRGetStandardErrorVec(KSP ksp,Vec *se)
283: {
284:   KSP_LSQR *lsqr = (KSP_LSQR*)ksp->data;

287:   *se = lsqr->se;
288:   return(0);
289: }

293: PetscErrorCode  KSPLSQRGetArnorm(KSP ksp,PetscReal *arnorm, PetscReal *rhs_norm, PetscReal *anorm)
294: {
295:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;

299:   *arnorm = lsqr->arnorm;
300:   if (anorm) {
301:     if (lsqr->anorm < 0.0) {
302:       PC  pc;
303:       Mat Amat;
304:       KSPGetPC(ksp,&pc);
305:       PCGetOperators(pc,&Amat,NULL);
306:       MatNorm(Amat,NORM_FROBENIUS,&lsqr->anorm);
307:     }
308:     *anorm = lsqr->anorm;
309:   }
310:   if (rhs_norm) *rhs_norm = lsqr->rhs_norm;
311:   return(0);
312: }

316: /*@C
317:    KSPLSQRMonitorDefault - Print the residual norm at each iteration of the LSQR method and the norm of the residual of the normal equations A'*A x = A' b

319:    Collective on KSP

321:    Input Parameters:
322: +  ksp   - iterative context
323: .  n     - iteration number
324: .  rnorm - 2-norm (preconditioned) residual value (may be estimated).
325: -  dummy - viewer and format context

327:    Level: intermediate

329: .keywords: KSP, default, monitor, residual

331: .seealso: KSPMonitorSet(), KSPMonitorTrueResidualNorm(), KSPMonitorLGResidualNormCreate(), KSPMonitorDefault()
332: @*/
333: PetscErrorCode  KSPLSQRMonitorDefault(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat *dummy)
334: {
336:   PetscViewer    viewer = dummy->viewer;
337:   KSP_LSQR       *lsqr  = (KSP_LSQR*)ksp->data;

340:   PetscViewerPushFormat(viewer,dummy->format);
341:   PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
342:   if (((PetscObject)ksp)->prefix) {
343:     PetscViewerASCIIPrintf(viewer,"  Residual norm and norm of normal equations for %s solve.\n",((PetscObject)ksp)->prefix);
344:   }
345:   if (!n) {
346:     PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm %14.12e\n",n,rnorm);
347:   } else {
348:     PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm %14.12e Residual norm normal equations %14.12e\n",n,rnorm,lsqr->arnorm);
349:   }
350:   PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
351:   PetscViewerPopFormat(viewer);
352:   return(0);
353: }

357: PetscErrorCode KSPSetFromOptions_LSQR(PetscOptionItems *PetscOptionsObject,KSP ksp)
358: {
360:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;

363:   PetscOptionsHead(PetscOptionsObject,"KSP LSQR Options");
364:   PetscOptionsName("-ksp_lsqr_set_standard_error","Set Standard Error Estimates of Solution","KSPLSQRSetStandardErrorVec",&lsqr->se_flg);
365:   KSPMonitorSetFromOptions(ksp,"-ksp_lsqr_monitor","Monitor residual norm and norm of residual of normal equations","KSPMonitorSet",KSPLSQRMonitorDefault);
366:   PetscOptionsTail();
367:   return(0);
368: }

372: PetscErrorCode KSPView_LSQR(KSP ksp,PetscViewer viewer)
373: {
374:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;
376:   PetscBool      iascii;

379:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
380:   if (iascii) {
381:     if (lsqr->se) {
382:       PetscReal rnorm;
383:       KSPLSQRGetStandardErrorVec(ksp,&lsqr->se);
384:       VecNorm(lsqr->se,NORM_2,&rnorm);
385:       PetscViewerASCIIPrintf(viewer,"  Norm of Standard Error %g, Iterations %D\n",(double)rnorm,ksp->its);
386:     }
387:   }
388:   return(0);
389: }

393: /*@C
394:    KSPLSQRDefaultConverged - Determines convergence of the LSQR Krylov method. This calls KSPConvergedDefault() and if that does not determine convergence then checks
395:       convergence for the least squares problem.

397:    Collective on KSP

399:    Input Parameters:
400: +  ksp   - iterative context
401: .  n     - iteration number
402: .  rnorm - 2-norm residual value (may be estimated)
403: -  ctx - convergence context which must be created by KSPConvergedDefaultCreate()

405:    reason is set to:
406: +   positive - if the iteration has converged;
407: .   negative - if residual norm exceeds divergence threshold;
408: -   0 - otherwise.

410:    Notes:
411:       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.

413:    Level: intermediate

415: .keywords: KSP, default, convergence, residual

417: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPConvergedSkip(), KSPConvergedReason, KSPGetConvergedReason(),
418:           KSPConvergedDefaultSetUIRNorm(), KSPConvergedDefaultSetUMIRNorm(), KSPConvergedDefaultCreate(), KSPConvergedDefaultDestroy(), KSPConvergedDefault()
419: @*/
420: PetscErrorCode  KSPLSQRDefaultConverged(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *ctx)
421: {
423:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;

426:   KSPConvergedDefault(ksp,n,rnorm,reason,ctx);
427:   if (!n || *reason) return(0);
428:   if (lsqr->arnorm/lsqr->rhs_norm < ksp->rtol) *reason = KSP_CONVERGED_RTOL_NORMAL;
429:   if (lsqr->arnorm < ksp->abstol) *reason = KSP_CONVERGED_ATOL_NORMAL;
430:   return(0);
431: }

433: /*MC
434:      KSPLSQR - This implements LSQR

436:    Options Database Keys:
437: +   -ksp_lsqr_set_standard_error  - Set Standard Error Estimates of Solution see KSPLSQRSetStandardErrorVec()
438: .   -ksp_lsqr_monitor - Monitor residual norm and norm of residual of normal equations
439: -   see KSPSolve()

441:    Level: beginner

443:    Notes:
444:      Supports non-square (rectangular) matrices.

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

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

452:      Supports only left preconditioning.

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

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

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


465:    For least squares problems without a zero to A*x = b, there are additional convergence tests for the residual of the normal equations, A'*(b - Ax), see KSPLSQRDefaultConverged()

467: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPLSQRDefaultConverged()

469: M*/
472: PETSC_EXTERN PetscErrorCode KSPCreate_LSQR(KSP ksp)
473: {
474:   KSP_LSQR       *lsqr;

478:   PetscNewLog(ksp,&lsqr);
479:   lsqr->se     = NULL;
480:   lsqr->se_flg = PETSC_FALSE;
481:   lsqr->arnorm = 0.0;
482:   ksp->data    = (void*)lsqr;
483:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,3);

485:   ksp->ops->setup          = KSPSetUp_LSQR;
486:   ksp->ops->solve          = KSPSolve_LSQR;
487:   ksp->ops->destroy        = KSPDestroy_LSQR;
488:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
489:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
490:   ksp->ops->setfromoptions = KSPSetFromOptions_LSQR;
491:   ksp->ops->view           = KSPView_LSQR;
492:   ksp->converged           = KSPLSQRDefaultConverged;
493:   return(0);
494: }