Actual source code: lsqr.c

  1: #define PETSCKSP_DLL

  3: /* lourens.vanzanen@shell.com contributed the standard error estimates of the solution, Jul 25, 2006 */

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

 7:  #include 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:   PetscTruth se_flg;   /* flag for -ksp_lsqr_set_standard_error */
 15: } KSP_LSQR;


 21: static PetscErrorCode KSPSetUp_LSQR(KSP ksp)
 22: {
 24:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;

 27:   if (ksp->pc_side == PC_SYMMETRIC){
 28:     SETERRQ(PETSC_ERR_SUP,"no symmetric preconditioning for KSPLSQR");
 29:   }

 31:   /* Get work vectors */
 32:   lsqr->nwork_m = 2;
 33:   if (lsqr->vwork_m) {
 34:     VecDestroyVecs(lsqr->vwork_m,lsqr->nwork_m);
 35:   }
 36:   lsqr->nwork_n = 4;
 37:   if (lsqr->vwork_n) {
 38:     VecDestroyVecs(lsqr->vwork_n,lsqr->nwork_n);
 39:   }
 40:   KSPGetVecs(ksp,lsqr->nwork_n,&lsqr->vwork_n,lsqr->nwork_m,&lsqr->vwork_m);
 41:   if (lsqr->se_flg && !lsqr->se){
 42:     /* lsqr->se is not set by user, get it from pmat */
 43:     Mat pmat;
 44:     PCGetOperators(ksp->pc,PETSC_NULL,&pmat,PETSC_NULL);
 45:     MatGetVecs(pmat,&lsqr->se,PETSC_NULL);
 46:   }
 47:   return(0);
 48: }

 52: static PetscErrorCode KSPSolve_LSQR(KSP ksp)
 53: {
 55:   PetscInt       i,size1,size2;
 56:   PetscScalar    rho,rhobar,phi,phibar,theta,c,s,tmp;
 57:   PetscReal      beta,alpha,rnorm;
 58:   Vec            X,B,V,V1,U,U1,TMP,W,W2,SE;
 59:   Mat            Amat,Pmat;
 60:   MatStructure   pflag;
 61:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;
 62:   PetscTruth     diagonalscale;

 65:   PCDiagonalScale(ksp->pc,&diagonalscale);
 66:   if (diagonalscale) SETERRQ1(PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);

 68:   PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);

 70:   /* vectors of length m, where system size is mxn */
 71:   B        = ksp->vec_rhs;
 72:   U        = lsqr->vwork_m[0];
 73:   U1       = lsqr->vwork_m[1];

 75:   /* vectors of length n */
 76:   X        = ksp->vec_sol;
 77:   W        = lsqr->vwork_n[0];
 78:   V        = lsqr->vwork_n[1];
 79:   V1       = lsqr->vwork_n[2];
 80:   W2       = lsqr->vwork_n[3];

 82:   /* standard error vector */
 83:   SE = lsqr->se;
 84:   if (SE){
 85:     VecGetSize(SE,&size1);
 86:     VecGetSize(X ,&size2);
 87:     if (size1 != size2) SETERRQ2(PETSC_ERR_ARG_SIZ,"Standard error vector (size %d) does not match solution vector (size %d)",size1,size2);
 88:     VecSet(SE,0.0);
 89:   }

 91:   /* Compute initial residual, temporarily use work vector u */
 92:   if (!ksp->guess_zero) {
 93:     KSP_MatMult(ksp,Amat,X,U);       /*   u <- b - Ax     */
 94:     VecAYPX(U,-1.0,B);
 95:   } else {
 96:     VecCopy(B,U);            /*   u <- b (x is 0) */
 97:   }

 99:   /* Test for nothing to do */
100:   VecNorm(U,NORM_2,&rnorm);
101:   PetscObjectTakeAccess(ksp);
102:   ksp->its   = 0;
103:   ksp->rnorm = rnorm;
104:   PetscObjectGrantAccess(ksp);
105:   KSPLogResidualHistory(ksp,rnorm);
106:   KSPMonitor(ksp,0,rnorm);
107:   (*ksp->converged)(ksp,0,rnorm,&ksp->reason,ksp->cnvP);
108:   if (ksp->reason) return(0);

110:   VecCopy(B,U);
111:   VecNorm(U,NORM_2,&beta);
112:   VecScale(U,1.0/beta);
113:   KSP_MatMultTranspose(ksp,Amat,U,V);
114:   VecNorm(V,NORM_2,&alpha);
115:   VecScale(V,1.0/alpha);

117:   VecCopy(V,W);
118:   VecSet(X,0.0);

120:   phibar = beta;
121:   rhobar = alpha;
122:   i = 0;
123:   do {

125:     KSP_MatMult(ksp,Amat,V,U1);
126:     VecAXPY(U1,-alpha,U);
127:     VecNorm(U1,NORM_2,&beta);
128:     VecScale(U1,1.0/beta);

130:     KSP_MatMultTranspose(ksp,Amat,U1,V1);
131:     VecAXPY(V1,-beta,V);
132:     VecNorm(V1,NORM_2,&alpha);
133:     VecScale(V1,1.0/alpha);

135:     rho    = PetscSqrtScalar(rhobar*rhobar + beta*beta);
136:     c      = rhobar / rho;
137:     s      = beta / rho;
138:     theta  = s * alpha;
139:     rhobar = - c * alpha;
140:     phi    = c * phibar;
141:     phibar = s * phibar;

143:     VecAXPY(X,phi/rho,W);  /*    x <- x + (phi/rho) w   */
144:     if (SE) {
145:       VecCopy(W,W2);
146:       VecSquare(W2);
147:       VecScale(W2,1.0/(rho*rho));
148:       VecAXPY(SE, 1.0, W2); /* SE <- SE + (w^2/rho^2) */
149:     }
150:     VecAYPX(W,-theta/rho,V1); /* w <- v - (theta/rho) w */

152:     rnorm = PetscRealPart(phibar);

154:     PetscObjectTakeAccess(ksp);
155:     ksp->its++;
156:     ksp->rnorm = rnorm;
157:     PetscObjectGrantAccess(ksp);
158:     KSPLogResidualHistory(ksp,rnorm);
159:     KSPMonitor(ksp,i+1,rnorm);
160:     (*ksp->converged)(ksp,i+1,rnorm,&ksp->reason,ksp->cnvP);
161:     if (ksp->reason) break;
162:     SWAP(U1,U,TMP);
163:     SWAP(V1,V,TMP);

165:     i++;
166:   } while (i<ksp->max_it);
167:   if (i >= ksp->max_it && !ksp->reason) {
168:     ksp->reason = KSP_DIVERGED_ITS;
169:   }

171:   /* Finish off the standard error estimates */
172:   if (SE) {
173:     tmp = 1.0;
174:     MatGetSize(Amat,&size1,&size2);
175:     if ( size1 > size2 ) tmp = size1 - size2;
176:     tmp = rnorm / PetscSqrtScalar(tmp);
177:     VecSqrt(SE);
178:     VecScale(SE,tmp);
179:   }
180:   return(0);
181: }

185: PetscErrorCode KSPDestroy_LSQR(KSP ksp)
186: {
187:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;


192:   /* Free work vectors */
193:   if (lsqr->vwork_n) {
194:     VecDestroyVecs(lsqr->vwork_n,lsqr->nwork_n);
195:   }
196:   if (lsqr->vwork_m) {
197:     VecDestroyVecs(lsqr->vwork_m,lsqr->nwork_m);
198:   }
199:   if (lsqr->se_flg && lsqr->se){
200:     VecDestroy(lsqr->se);
201:   }
202:   PetscFree(ksp->data);
203:   return(0);
204: }

208: PetscErrorCode  KSPLSQRSetStandardErrorVec( KSP ksp, Vec se )
209: {
210:   KSP_LSQR  *lsqr = (KSP_LSQR*)ksp->data;

214:   if (lsqr->se) {
215:     VecDestroy(lsqr->se);
216:   }
217:   lsqr->se     = se;
218:   return(0);
219: }

223: PetscErrorCode  KSPLSQRGetStandardErrorVec( KSP ksp,Vec *se )
224: {
225:   KSP_LSQR *lsqr = (KSP_LSQR*)ksp->data;

228:   *se = lsqr->se;
229:   return(0);
230: }

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

240:   PetscOptionsHead("KSP LSQR Options");
241:   PetscOptionsName("-ksp_lsqr_set_standard_error","Set Standard Error Estimates of Solution","KSPLSQRSetStandardErrorVec",&lsqr->se_flg);
242:   PetscOptionsTail();
243:   return(0);
244: }

248: PetscErrorCode KSPView_LSQR(KSP ksp,PetscViewer viewer)
249: {
250:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;

254:   if (lsqr->se) {
255:     PetscReal rnorm;
256:     KSPLSQRGetStandardErrorVec(ksp,&lsqr->se);
257:     VecNorm(lsqr->se,NORM_2,&rnorm);
258:     PetscPrintf(PETSC_COMM_WORLD,"  Norm of Standard Error %A, Iterations %D\n",rnorm,ksp->its);
259:   }
260:   return(0);
261: }

263: /*MC
264:      KSPLSQR - This implements LSQR

266:    Options Database Keys:
267: .   see KSPSolve()

269:    Level: beginner

271:    Notes:  This algorithm DOES NOT use a preconditioner. It ignores any preconditioner arguments specified.
272:            Reference: Paige and Saunders, ACM Transactions on Mathematical Software, Vol 8, pp 43-71, 1982

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

276: M*/
280: PetscErrorCode  KSPCreate_LSQR(KSP ksp)
281: {
282:   KSP_LSQR       *lsqr;

286:   PetscNewLog(ksp,KSP_LSQR,&lsqr);
287:   lsqr->se     = PETSC_NULL;
288:   lsqr->se_flg = PETSC_FALSE;
289:   PCSetType(ksp->pc,PCNONE);
290:   ksp->data                      = (void*)lsqr;
291:   ksp->pc_side                   = PC_LEFT;
292:   ksp->ops->setup                = KSPSetUp_LSQR;
293:   ksp->ops->solve                = KSPSolve_LSQR;
294:   ksp->ops->destroy              = KSPDestroy_LSQR;
295:   ksp->ops->buildsolution        = KSPDefaultBuildSolution;
296:   ksp->ops->buildresidual        = KSPDefaultBuildResidual;
297:   ksp->ops->setfromoptions       = KSPSetFromOptions_LSQR;
298:   ksp->ops->view                 = KSPView_LSQR;
299:   return(0);
300: }

305: PetscErrorCode  VecSquare(Vec v)
306: {
308:   PetscScalar    *x;
309:   PetscInt       i, n;

312:   VecGetLocalSize(v, &n);
313:   VecGetArray(v, &x);
314:   for(i = 0; i < n; i++) {
315:     x[i] *= x[i];
316:   }
317:   VecRestoreArray(v, &x);
318:   return(0);
319: }