Actual source code: lsqr.c

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

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

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

 14: static PetscErrorCode KSPSetUp_LSQR(KSP ksp)
 15: {
 17:   PetscInt       nw;
 18:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;

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

 25:   /* Get work vectors */
 26:   lsqr->nwork_m = nw = 2;
 27:   if (lsqr->vwork_m) {
 28:     VecDestroyVecs(lsqr->vwork_m,lsqr->nwork_m);
 29:   }
 30:   KSPGetVecs(ksp,nw,&lsqr->vwork_m);
 31:   PetscLogObjectParents(ksp,nw,lsqr->vwork_m);

 33:   lsqr->nwork_n = nw = 3;
 34:   if (lsqr->vwork_n) {
 35:     VecDestroyVecs(lsqr->vwork_n,lsqr->nwork_n);
 36:   }
 37:   KSPGetVecs(ksp,nw,&lsqr->vwork_n);
 38:   PetscLogObjectParents(ksp,nw,lsqr->vwork_n);

 40:   return(0);
 41: }

 45: static PetscErrorCode KSPSolve_LSQR(KSP ksp)
 46: {
 48:   PetscInt       i;
 49:   PetscScalar    rho,rhobar,phi,phibar,theta,c,s,tmp,zero = 0.0,mone=-1.0;
 50:   PetscReal      beta,alpha,rnorm;
 51:   Vec            X,B,V,V1,U,U1,TMP,W;
 52:   Mat            Amat,Pmat;
 53:   MatStructure   pflag;
 54:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;
 55:   PetscTruth     diagonalscale;

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

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

 63:   /* vectors of length m, where system size is mxn */
 64:   B        = ksp->vec_rhs;
 65:   U        = lsqr->vwork_m[0];
 66:   U1       = lsqr->vwork_m[1];

 68:   /* vectors of length n */
 69:   X        = ksp->vec_sol;
 70:   W        = lsqr->vwork_n[0];
 71:   V        = lsqr->vwork_n[1];
 72:   V1       = lsqr->vwork_n[2];

 74:   /* Compute initial residual, temporarily use work vector u */
 75:   if (!ksp->guess_zero) {
 76:     KSP_MatMult(ksp,Amat,X,U);       /*   u <- b - Ax     */
 77:     VecAYPX(&mone,B,U);
 78:   } else {
 79:     VecCopy(B,U);            /*   u <- b (x is 0) */
 80:   }

 82:   /* Test for nothing to do */
 83:   VecNorm(U,NORM_2,&rnorm);
 84:   PetscObjectTakeAccess(ksp);
 85:   ksp->its   = 0;
 86:   ksp->rnorm = rnorm;
 87:   PetscObjectGrantAccess(ksp);
 88:   KSPLogResidualHistory(ksp,rnorm);
 89:   KSPMonitor(ksp,0,rnorm);
 90:   (*ksp->converged)(ksp,0,rnorm,&ksp->reason,ksp->cnvP);
 91:   if (ksp->reason) return(0);

 93:   VecCopy(B,U);
 94:   VecNorm(U,NORM_2,&beta);
 95:   tmp = 1.0/beta; VecScale(&tmp,U);
 96:   KSP_MatMultTranspose(ksp,Amat,U,V);
 97:   VecNorm(V,NORM_2,&alpha);
 98:   tmp = 1.0/alpha; VecScale(&tmp,V);

100:   VecCopy(V,W);
101:   VecSet(&zero,X);

103:   phibar = beta;
104:   rhobar = alpha;
105:   i = 0;
106:   do {

108:     KSP_MatMult(ksp,Amat,V,U1);
109:     tmp  = -alpha; VecAXPY(&tmp,U,U1);
110:     VecNorm(U1,NORM_2,&beta);
111:     tmp  = 1.0/beta; VecScale(&tmp,U1);

113:     KSP_MatMultTranspose(ksp,Amat,U1,V1);
114:     tmp  = -beta; VecAXPY(&tmp,V,V1);
115:     VecNorm(V1,NORM_2,&alpha);
116:     tmp  = 1.0 / alpha; VecScale(&tmp,V1);

118:     rho    = PetscSqrtScalar(rhobar*rhobar + beta*beta);
119:     c      = rhobar / rho;
120:     s      = beta / rho;
121:     theta  = s * alpha;
122:     rhobar = - c * alpha;
123:     phi    = c * phibar;
124:     phibar = s * phibar;

126:     tmp  = phi/rho;
127:     VecAXPY(&tmp,W,X);  /*    x <- x + (phi/rho) w   */
128:     tmp  = -theta/rho;
129:     VecAYPX(&tmp,V1,W); /*    w <- v - (theta/rho) w */

131:     rnorm = PetscRealPart(phibar);

133:     PetscObjectTakeAccess(ksp);
134:     ksp->its++;
135:     ksp->rnorm = rnorm;
136:     PetscObjectGrantAccess(ksp);
137:     KSPLogResidualHistory(ksp,rnorm);
138:     KSPMonitor(ksp,i+1,rnorm);
139:     (*ksp->converged)(ksp,i+1,rnorm,&ksp->reason,ksp->cnvP);
140:     if (ksp->reason) break;
141:     SWAP(U1,U,TMP);
142:     SWAP(V1,V,TMP);

144:     i++;
145:   } while (i<ksp->max_it);
146:   if (i >= ksp->max_it && !ksp->reason) {
147:     ksp->reason = KSP_DIVERGED_ITS;
148:   }

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

152:   return(0);
153: }

157: PetscErrorCode KSPDestroy_LSQR(KSP ksp)
158: {
159:   KSP_LSQR       *lsqr = (KSP_LSQR*)ksp->data;


164:   /* Free work vectors */
165:   if (lsqr->vwork_n) {
166:     VecDestroyVecs(lsqr->vwork_n,lsqr->nwork_n);
167:   }
168:   if (lsqr->vwork_m) {
169:     VecDestroyVecs(lsqr->vwork_m,lsqr->nwork_m);
170:   }
171:   PetscFree(lsqr);
172:   return(0);
173: }

175: /*MC
176:      KSPLSQR - This implements LSQR

178:    Options Database Keys:
179: .   see KSPSolve()

181:    Level: beginner

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

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

188: M*/
192: PetscErrorCode KSPCreate_LSQR(KSP ksp)
193: {
194:   KSP_LSQR       *lsqr;

198:   PetscMalloc(sizeof(KSP_LSQR),&lsqr);
199:   PetscMemzero(lsqr,sizeof(KSP_LSQR));
200:   PetscLogObjectMemory(ksp,sizeof(KSP_LSQR));
201:   ksp->data                      = (void*)lsqr;
202:   ksp->pc_side                   = PC_LEFT;
203:   ksp->ops->setup                = KSPSetUp_LSQR;
204:   ksp->ops->solve                = KSPSolve_LSQR;
205:   ksp->ops->destroy              = KSPDestroy_LSQR;
206:   ksp->ops->buildsolution        = KSPDefaultBuildSolution;
207:   ksp->ops->buildresidual        = KSPDefaultBuildResidual;
208:   ksp->ops->setfromoptions       = 0;
209:   ksp->ops->view                 = 0;
210:   return(0);
211: }