Actual source code: lsqr.c

  1: #define PETSCKSP_DLL

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

 5:  #include src/ksp/ksp/kspimpl.h

  7: typedef struct {
  8:   PetscInt  nwork_n,nwork_m;
  9:   Vec       *vwork_m;  /* work vectors of length m, where the system is size m x n */
 10:   Vec       *vwork_n;  /* work vectors of length n */
 11: } KSP_LSQR;

 15: static PetscErrorCode KSPSetUp_LSQR(KSP ksp)
 16: {
 18:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;

 21:   if (ksp->pc_side == PC_SYMMETRIC){
 22:     SETERRQ(PETSC_ERR_SUP,"no symmetric preconditioning for KSPLSQR");
 23:   }

 25:   /* Get work vectors */
 26:   lsqr->nwork_m = 2;
 27:   if (lsqr->vwork_m) {
 28:     VecDestroyVecs(lsqr->vwork_m,lsqr->nwork_m);
 29:   }
 30:   lsqr->nwork_n = 3;
 31:   if (lsqr->vwork_n) {
 32:     VecDestroyVecs(lsqr->vwork_n,lsqr->nwork_n);
 33:   }
 34:   KSPGetVecs(ksp,lsqr->nwork_n,&lsqr->vwork_n,lsqr->nwork_m,&lsqr->vwork_m);
 35:   return(0);
 36: }

 40: static PetscErrorCode KSPSolve_LSQR(KSP ksp)
 41: {
 43:   PetscInt       i;
 44:   PetscScalar    rho,rhobar,phi,phibar,theta,c,s;
 45:   PetscReal      beta,alpha,rnorm;
 46:   Vec            X,B,V,V1,U,U1,TMP,W;
 47:   Mat            Amat,Pmat;
 48:   MatStructure   pflag;
 49:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;
 50:   PetscTruth     diagonalscale;

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

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

 58:   /* vectors of length m, where system size is mxn */
 59:   B        = ksp->vec_rhs;
 60:   U        = lsqr->vwork_m[0];
 61:   U1       = lsqr->vwork_m[1];

 63:   /* vectors of length n */
 64:   X        = ksp->vec_sol;
 65:   W        = lsqr->vwork_n[0];
 66:   V        = lsqr->vwork_n[1];
 67:   V1       = lsqr->vwork_n[2];

 69:   /* Compute initial residual, temporarily use work vector u */
 70:   if (!ksp->guess_zero) {
 71:     KSP_MatMult(ksp,Amat,X,U);       /*   u <- b - Ax     */
 72:     VecAYPX(U,-1.0,B);
 73:   } else {
 74:     VecCopy(B,U);            /*   u <- b (x is 0) */
 75:   }

 77:   /* Test for nothing to do */
 78:   VecNorm(U,NORM_2,&rnorm);
 79:   PetscObjectTakeAccess(ksp);
 80:   ksp->its   = 0;
 81:   ksp->rnorm = rnorm;
 82:   PetscObjectGrantAccess(ksp);
 83:   KSPLogResidualHistory(ksp,rnorm);
 84:   KSPMonitor(ksp,0,rnorm);
 85:   (*ksp->converged)(ksp,0,rnorm,&ksp->reason,ksp->cnvP);
 86:   if (ksp->reason) return(0);

 88:   VecCopy(B,U);
 89:   VecNorm(U,NORM_2,&beta);
 90:   VecScale(U,1.0/beta);
 91:   KSP_MatMultTranspose(ksp,Amat,U,V);
 92:   VecNorm(V,NORM_2,&alpha);
 93:   VecScale(V,1.0/alpha);

 95:   VecCopy(V,W);
 96:   VecSet(X,0.0);

 98:   phibar = beta;
 99:   rhobar = alpha;
100:   i = 0;
101:   do {

103:     KSP_MatMult(ksp,Amat,V,U1);
104:     VecAXPY(U1,-alpha,U);
105:     VecNorm(U1,NORM_2,&beta);
106:     VecScale(U1,1.0/beta);

108:     KSP_MatMultTranspose(ksp,Amat,U1,V1);
109:     VecAXPY(V1,-beta,V);
110:     VecNorm(V1,NORM_2,&alpha);
111:     VecScale(V1,1.0/alpha);

113:     rho    = PetscSqrtScalar(rhobar*rhobar + beta*beta);
114:     c      = rhobar / rho;
115:     s      = beta / rho;
116:     theta  = s * alpha;
117:     rhobar = - c * alpha;
118:     phi    = c * phibar;
119:     phibar = s * phibar;

121:     VecAXPY(X,phi/rho,W);  /*    x <- x + (phi/rho) w   */
122:     VecAYPX(W,-theta/rho,V1); /*    w <- v - (theta/rho) w */

124:     rnorm = PetscRealPart(phibar);

126:     PetscObjectTakeAccess(ksp);
127:     ksp->its++;
128:     ksp->rnorm = rnorm;
129:     PetscObjectGrantAccess(ksp);
130:     KSPLogResidualHistory(ksp,rnorm);
131:     KSPMonitor(ksp,i+1,rnorm);
132:     (*ksp->converged)(ksp,i+1,rnorm,&ksp->reason,ksp->cnvP);
133:     if (ksp->reason) break;
134:     SWAP(U1,U,TMP);
135:     SWAP(V1,V,TMP);

137:     i++;
138:   } while (i<ksp->max_it);
139:   if (i >= ksp->max_it && !ksp->reason) {
140:     ksp->reason = KSP_DIVERGED_ITS;
141:   }

143:   /* KSPUnwindPreconditioner(ksp,X,W); */

145:   return(0);
146: }

150: PetscErrorCode KSPDestroy_LSQR(KSP ksp)
151: {
152:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;


157:   /* Free work vectors */
158:   if (lsqr->vwork_n) {
159:     VecDestroyVecs(lsqr->vwork_n,lsqr->nwork_n);
160:   }
161:   if (lsqr->vwork_m) {
162:     VecDestroyVecs(lsqr->vwork_m,lsqr->nwork_m);
163:   }
164:   PetscFree(lsqr);
165:   return(0);
166: }

168: /*MC
169:      KSPLSQR - This implements LSQR

171:    Options Database Keys:
172: .   see KSPSolve()

174:    Level: beginner

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

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

181: M*/
185: PetscErrorCode PETSCKSP_DLLEXPORT KSPCreate_LSQR(KSP ksp)
186: {
187:   KSP_LSQR       *lsqr;

191:   PetscNew(KSP_LSQR,&lsqr);
192:   PCSetType(ksp->pc,PCNONE);
193:   ksp->data                      = (void*)lsqr;
194:   ksp->pc_side                   = PC_LEFT;
195:   ksp->ops->setup                = KSPSetUp_LSQR;
196:   ksp->ops->solve                = KSPSolve_LSQR;
197:   ksp->ops->destroy              = KSPDestroy_LSQR;
198:   ksp->ops->buildsolution        = KSPDefaultBuildSolution;
199:   ksp->ops->buildresidual        = KSPDefaultBuildResidual;
200:   ksp->ops->setfromoptions       = 0;
201:   ksp->ops->view                 = 0;
202:   return(0);
203: }