Actual source code: lsqr.c

petsc-3.6.1 2015-08-06
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>
  8: #include <../src/ksp/ksp/impls/lsqr/lsqr.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:   PetscReal arnorm;     /* Norm of the vector A.r */
 17:   PetscReal anorm;      /* Frobenius norm of the matrix A */
 18:   PetscReal rhs_norm;   /* Norm of the right hand side */
 19: } KSP_LSQR;

 21: extern PetscErrorCode  VecSquare(Vec);

 25: static PetscErrorCode KSPSetUp_LSQR(KSP ksp)
 26: {
 28:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;
 29:   PetscBool      nopreconditioner;

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

 35:   lsqr->nwork_m = 2;
 36:   if (lsqr->vwork_m) {
 37:     VecDestroyVecs(lsqr->nwork_m,&lsqr->vwork_m);
 38:   }
 39:   if (nopreconditioner) lsqr->nwork_n = 4;
 40:   else lsqr->nwork_n = 5;

 42:   if (lsqr->vwork_n) {
 43:     VecDestroyVecs(lsqr->nwork_n,&lsqr->vwork_n);
 44:   }
 45:   KSPCreateVecs(ksp,lsqr->nwork_n,&lsqr->vwork_n,lsqr->nwork_m,&lsqr->vwork_m);
 46:   if (lsqr->se_flg && !lsqr->se) {
 47:     /* lsqr->se is not set by user, get it from pmat */
 48:     Vec *se;
 49:     KSPCreateVecs(ksp,1,&se,0,NULL);
 50:     lsqr->se = *se;
 51:     PetscFree(se);
 52:   }
 53:   return(0);
 54: }

 58: static PetscErrorCode KSPSolve_LSQR(KSP ksp)
 59: {
 61:   PetscInt       i,size1,size2;
 62:   PetscScalar    rho,rhobar,phi,phibar,theta,c,s,tmp,tau;
 63:   PetscReal      beta,alpha,rnorm;
 64:   Vec            X,B,V,V1,U,U1,TMP,W,W2,SE,Z = NULL;
 65:   Mat            Amat,Pmat;
 66:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;
 67:   PetscBool      diagonalscale,nopreconditioner;

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

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

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

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

 83:   /* vectors of length m, where system size is mxn */
 84:   B  = ksp->vec_rhs;
 85:   U  = lsqr->vwork_m[0];
 86:   U1 = lsqr->vwork_m[1];

 88:   /* vectors of length n */
 89:   X  = ksp->vec_sol;
 90:   W  = lsqr->vwork_n[0];
 91:   V  = lsqr->vwork_n[1];
 92:   V1 = lsqr->vwork_n[2];
 93:   W2 = lsqr->vwork_n[3];
 94:   if (!nopreconditioner) Z = lsqr->vwork_n[4];

 96:   /* standard error vector */
 97:   SE = lsqr->se;
 98:   if (SE) {
 99:     VecGetSize(SE,&size1);
100:     VecGetSize(X,&size2);
101:     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);
102:     VecSet(SE,0.0);
103:   }

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

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

124:   beta = rnorm;
125:   VecScale(U,1.0/beta);
126:   KSP_MatMultTranspose(ksp,Amat,U,V);
127:   if (nopreconditioner) {
128:     VecNorm(V,NORM_2,&alpha);
129:   } else {
130:     PCApply(ksp->pc,V,Z);
131:     VecDotRealPart(V,Z,&alpha);
132:     if (alpha <= 0.0) {
133:       ksp->reason = KSP_DIVERGED_BREAKDOWN;
134:       return(0);
135:     }
136:     alpha = PetscSqrtReal(alpha);
137:     VecScale(Z,1.0/alpha);
138:   }
139:   VecScale(V,1.0/alpha);

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

147:   lsqr->arnorm = alpha * beta;
148:   phibar       = beta;
149:   rhobar       = alpha;
150:   i            = 0;
151:   do {
152:     if (nopreconditioner) {
153:       KSP_MatMult(ksp,Amat,V,U1);
154:     } else {
155:       KSP_MatMult(ksp,Amat,Z,U1);
156:     }
157:     VecAXPY(U1,-alpha,U);
158:     VecNorm(U1,NORM_2,&beta);
159:     if (beta > 0.0) {
160:       VecScale(U1,1.0/beta); /* beta*U1 = Amat*V - alpha*U */
161:     }

163:     KSP_MatMultTranspose(ksp,Amat,U1,V1);
164:     VecAXPY(V1,-beta,V);
165:     if (nopreconditioner) {
166:       VecNorm(V1,NORM_2,&alpha);
167:     } else {
168:       PCApply(ksp->pc,V1,Z);
169:       VecDotRealPart(V1,Z,&alpha);
170:       if (alpha <= 0.0) {
171:         ksp->reason = KSP_DIVERGED_BREAKDOWN;
172:         break;
173:       }
174:       alpha = PetscSqrtReal(alpha);
175:       VecScale(Z,1.0/alpha);
176:     }
177:     VecScale(V1,1.0/alpha); /* alpha*V1 = Amat^T*U1 - beta*V */
178:     rho    = PetscSqrtScalar(rhobar*rhobar + beta*beta);
179:     c      = rhobar / rho;
180:     s      = beta / rho;
181:     theta  = s * alpha;
182:     rhobar = -c * alpha;
183:     phi    = c * phibar;
184:     phibar = s * phibar;
185:     tau    = s * phi;

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

189:     if (SE) {
190:       VecCopy(W,W2);
191:       VecSquare(W2);
192:       VecScale(W2,1.0/(rho*rho));
193:       VecAXPY(SE, 1.0, W2); /* SE <- SE + (w^2/rho^2) */
194:     }
195:     if (nopreconditioner) {
196:       VecAYPX(W,-theta/rho,V1);  /* w <- v - (theta/rho) w */
197:     } else {
198:       VecAYPX(W,-theta/rho,Z);   /* w <- z - (theta/rho) w */
199:     }

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

204:     PetscObjectSAWsTakeAccess((PetscObject)ksp);
205:     ksp->its++;
206:     ksp->rnorm = rnorm;
207:     PetscObjectSAWsGrantAccess((PetscObject)ksp);
208:     KSPLogResidualHistory(ksp,rnorm);
209:     KSPMonitor(ksp,i+1,rnorm);
210:     (*ksp->converged)(ksp,i+1,rnorm,&ksp->reason,ksp->cnvP);
211:     if (ksp->reason) break;
212:     SWAP(U1,U,TMP);
213:     SWAP(V1,V,TMP);

215:     i++;
216:   } while (i<ksp->max_it);
217:   if (i >= ksp->max_it && !ksp->reason) ksp->reason = KSP_DIVERGED_ITS;

219:   /* Finish off the standard error estimates */
220:   if (SE) {
221:     tmp  = 1.0;
222:     MatGetSize(Amat,&size1,&size2);
223:     if (size1 > size2) tmp = size1 - size2;
224:     tmp  = rnorm / PetscSqrtScalar(tmp);
225:     VecSqrtAbs(SE);
226:     VecScale(SE,tmp);
227:   }
228:   return(0);
229: }


234: PetscErrorCode KSPDestroy_LSQR(KSP ksp)
235: {
236:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;

240:   /* Free work vectors */
241:   if (lsqr->vwork_n) {
242:     VecDestroyVecs(lsqr->nwork_n,&lsqr->vwork_n);
243:   }
244:   if (lsqr->vwork_m) {
245:     VecDestroyVecs(lsqr->nwork_m,&lsqr->vwork_m);
246:   }
247:   if (lsqr->se_flg) {
248:     VecDestroy(&lsqr->se);
249:   }
250:   PetscFree(ksp->data);
251:   return(0);
252: }

256: PetscErrorCode  KSPLSQRSetStandardErrorVec(KSP ksp, Vec se)
257: {
258:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;

262:   VecDestroy(&lsqr->se);
263:   lsqr->se = se;
264:   return(0);
265: }

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

274:   *se = lsqr->se;
275:   return(0);
276: }

280: PetscErrorCode  KSPLSQRGetArnorm(KSP ksp,PetscReal *arnorm, PetscReal *rhs_norm, PetscReal *anorm)
281: {
282:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;

286:   *arnorm = lsqr->arnorm;
287:   if (anorm) {
288:     if (lsqr->anorm < 0.0) {
289:       PC  pc;
290:       Mat Amat;
291:       KSPGetPC(ksp,&pc);
292:       PCGetOperators(pc,&Amat,NULL);
293:       MatNorm(Amat,NORM_FROBENIUS,&lsqr->anorm);
294:     }
295:     *anorm = lsqr->anorm;
296:   }
297:   if (rhs_norm) *rhs_norm = lsqr->rhs_norm;
298:   return(0);
299: }

303: /*@C
304:    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

306:    Collective on KSP

308:    Input Parameters:
309: +  ksp   - iterative context
310: .  n     - iteration number
311: .  rnorm - 2-norm (preconditioned) residual value (may be estimated).
312: -  dummy - unused monitor context

314:    Level: intermediate

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

318: .seealso: KSPMonitorSet(), KSPMonitorTrueResidualNorm(), KSPMonitorLGResidualNormCreate(), KSPMonitorDefault()
319: @*/
320: PetscErrorCode  KSPLSQRMonitorDefault(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy)
321: {
323:   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)ksp));
324:   KSP_LSQR       *lsqr  = (KSP_LSQR*)ksp->data;

327:   PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
328:   if (((PetscObject)ksp)->prefix) {
329:     PetscViewerASCIIPrintf(viewer,"  Residual norm and norm of normal equations for %s solve.\n",((PetscObject)ksp)->prefix);
330:   }
331:   if (!n) {
332:     PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm %14.12e\n",n,rnorm);
333:   } else {
334:     PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm %14.12e Residual norm normal equations %14.12e\n",n,rnorm,lsqr->arnorm);
335:   }
336:   PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
337:   return(0);
338: }

342: PetscErrorCode KSPSetFromOptions_LSQR(PetscOptions *PetscOptionsObject,KSP ksp)
343: {
345:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;
346:   char           monfilename[PETSC_MAX_PATH_LEN];
347:   PetscViewer    monviewer;
348:   PetscBool      flg;

351:   PetscOptionsHead(PetscOptionsObject,"KSP LSQR Options");
352:   PetscOptionsName("-ksp_lsqr_set_standard_error","Set Standard Error Estimates of Solution","KSPLSQRSetStandardErrorVec",&lsqr->se_flg);
353:   PetscOptionsString("-ksp_lsqr_monitor","Monitor residual norm and norm of residual of normal equations","KSPMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
354:   if (flg) {
355:     PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ksp),monfilename,&monviewer);
356:     KSPMonitorSet(ksp,KSPLSQRMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
357:   }
358:   PetscOptionsTail();
359:   return(0);
360: }

364: PetscErrorCode KSPView_LSQR(KSP ksp,PetscViewer viewer)
365: {
366:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;
368:   PetscBool      iascii;

371:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
372:   if (iascii) {
373:     if (lsqr->se) {
374:       PetscReal rnorm;
375:       KSPLSQRGetStandardErrorVec(ksp,&lsqr->se);
376:       VecNorm(lsqr->se,NORM_2,&rnorm);
377:       PetscViewerASCIIPrintf(viewer,"  Norm of Standard Error %g, Iterations %D\n",(double)rnorm,ksp->its);
378:     }
379:   }
380:   return(0);
381: }

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

389:    Collective on KSP

391:    Input Parameters:
392: +  ksp   - iterative context
393: .  n     - iteration number
394: .  rnorm - 2-norm residual value (may be estimated)
395: -  ctx - convergence context which must be created by KSPConvergedDefaultCreate()

397:    reason is set to:
398: +   positive - if the iteration has converged;
399: .   negative - if residual norm exceeds divergence threshold;
400: -   0 - otherwise.

402:    Notes:
403:       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.

405:    Level: intermediate

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

409: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPConvergedSkip(), KSPConvergedReason, KSPGetConvergedReason(),
410:           KSPConvergedDefaultSetUIRNorm(), KSPConvergedDefaultSetUMIRNorm(), KSPConvergedDefaultCreate(), KSPConvergedDefaultDestroy(), KSPConvergedDefault()
411: @*/
412: PetscErrorCode  KSPLSQRDefaultConverged(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *ctx)
413: {
415:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;

418:   KSPConvergedDefault(ksp,n,rnorm,reason,ctx);
419:   if (!n || *reason) return(0);
420:   if (lsqr->arnorm/lsqr->rhs_norm < ksp->rtol) *reason = KSP_CONVERGED_RTOL_NORMAL;
421:   if (lsqr->arnorm < ksp->abstol) *reason = KSP_CONVERGED_ATOL_NORMAL;
422:   return(0);
423: }



427: /*MC
428:      KSPLSQR - This implements LSQR

430:    Options Database Keys:
431: +   -ksp_lsqr_set_standard_error  - Set Standard Error Estimates of Solution see KSPLSQRSetStandardErrorVec()
432: .   -ksp_lsqr_monitor - Monitor residual norm and norm of residual of normal equations
433: -   see KSPSolve()

435:    Level: beginner

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

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

444:      Supports only left preconditioning.

446:    References:The original unpreconditioned algorithm can be found in Paige and Saunders, ACM Transactions on Mathematical Software, Vol 8, pp 43-71, 1982.
447:      In exact arithmetic the LSQR method (with no preconditioning) is identical to the KSPCG algorithm applied to the normal equations.
448:      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
449:      track the true norm of the residual and uses that in the convergence test.

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


455:    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()

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

459: M*/
462: PETSC_EXTERN PetscErrorCode KSPCreate_LSQR(KSP ksp)
463: {
464:   KSP_LSQR       *lsqr;

468:   PetscNewLog(ksp,&lsqr);
469:   lsqr->se     = NULL;
470:   lsqr->se_flg = PETSC_FALSE;
471:   lsqr->arnorm = 0.0;
472:   ksp->data    = (void*)lsqr;
473:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,3);

475:   ksp->ops->setup          = KSPSetUp_LSQR;
476:   ksp->ops->solve          = KSPSolve_LSQR;
477:   ksp->ops->destroy        = KSPDestroy_LSQR;
478:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
479:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
480:   ksp->ops->setfromoptions = KSPSetFromOptions_LSQR;
481:   ksp->ops->view           = KSPView_LSQR;
482:   ksp->converged           = KSPLSQRDefaultConverged;
483:   return(0);
484: }

488: PetscErrorCode  VecSquare(Vec v)
489: {
491:   PetscScalar    *x;
492:   PetscInt       i, n;

495:   VecGetLocalSize(v, &n);
496:   VecGetArray(v, &x);
497:   for (i = 0; i < n; i++) x[i] *= PetscConj(x[i]);
498:   VecRestoreArray(v, &x);
499:   return(0);
500: }