Actual source code: lsqr.c
1: /*$Id: lsqr.c,v 1.69 2001/08/07 03:03:53 balay Exp $*/
3: #define SWAP(a,b,c) { c = a; a = b; b = c; }
5: /*
6: This implements LSQR (Paige and Saunders, ACM Transactions on
7: Mathematical Software, Vol 8, pp 43-71, 1982).
9: This algorithm DOES NOT use a preconditioner. It ignores
10: any preconditioner arguments specified.
11: */
12: #include src/sles/ksp/kspimpl.h
14: typedef struct {
15: int nwork_n,nwork_m;
16: Vec *vwork_m; /* work vectors of length m, where the system is size m x n */
17: Vec *vwork_n; /* work vectors of length m */
18: } KSP_LSQR;
22: static int KSPSetUp_LSQR(KSP ksp)
23: {
24: int ierr,nw;
25: KSP_LSQR *lsqr = (KSP_LSQR*)ksp->data;
28: if (ksp->pc_side == PC_SYMMETRIC){
29: SETERRQ(2,"no symmetric preconditioning for KSPLSQR");
30: }
32: /* Get work vectors */
33: lsqr->nwork_m = nw = 2;
34: if (lsqr->vwork_m) {
35: VecDestroyVecs(lsqr->vwork_m,lsqr->nwork_m);
36: }
37: VecDuplicateVecs(ksp->vec_rhs,nw,&lsqr->vwork_m);
38: PetscLogObjectParents(ksp,nw,lsqr->vwork_m);
40: lsqr->nwork_n = nw = 3;
41: if (lsqr->vwork_n) {
42: VecDestroyVecs(lsqr->vwork_n,lsqr->nwork_n);
43: }
44: VecDuplicateVecs(ksp->vec_sol,nw,&lsqr->vwork_n);
45: PetscLogObjectParents(ksp,nw,lsqr->vwork_n);
47: return(0);
48: }
52: static int KSPSolve_LSQR(KSP ksp,int *its)
53: {
54: int i,maxit,ierr;
55: PetscScalar rho,rhobar,phi,phibar,theta,c,s,tmp,zero = 0.0,mone=-1.0;
56: PetscReal beta,alpha,rnorm;
57: Vec X,B,V,V1,U,U1,TMP,W;
58: Mat Amat,Pmat;
59: MatStructure pflag;
60: KSP_LSQR *lsqr = (KSP_LSQR*)ksp->data;
61: PetscTruth diagonalscale;
64: PCDiagonalScale(ksp->B,&diagonalscale);
65: if (diagonalscale) SETERRQ1(1,"Krylov method %s does not support diagonal scaling",ksp->type_name);
67: PCGetOperators(ksp->B,&Amat,&Pmat,&pflag);
68: maxit = ksp->max_it;
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];
81: /* Compute initial residual, temporarily use work vector u */
82: if (!ksp->guess_zero) {
83: KSP_MatMult(ksp,Amat,X,U); /* u <- b - Ax */
84: VecAYPX(&mone,B,U);
85: } else {
86: VecCopy(B,U); /* u <- b (x is 0) */
87: }
89: /* Test for nothing to do */
90: VecNorm(U,NORM_2,&rnorm);
91: PetscObjectTakeAccess(ksp);
92: ksp->its = 0;
93: ksp->rnorm = rnorm;
94: PetscObjectGrantAccess(ksp);
95: (*ksp->converged)(ksp,0,rnorm,&ksp->reason,ksp->cnvP);
96: if (ksp->reason) {*its = 0; return(0);}
97: KSPLogResidualHistory(ksp,rnorm);
98: KSPMonitor(ksp,0,rnorm);
100: VecCopy(B,U);
101: VecNorm(U,NORM_2,&beta);
102: tmp = 1.0/beta; VecScale(&tmp,U);
103: KSP_MatMultTranspose(ksp,Amat,U,V);
104: VecNorm(V,NORM_2,&alpha);
105: tmp = 1.0/alpha; VecScale(&tmp,V);
107: VecCopy(V,W);
108: VecSet(&zero,X);
110: phibar = beta;
111: rhobar = alpha;
112: for (i=0; i<maxit; i++) {
114: KSP_MatMult(ksp,Amat,V,U1);
115: tmp = -alpha; VecAXPY(&tmp,U,U1);
116: VecNorm(U1,NORM_2,&beta);
117: tmp = 1.0/beta; VecScale(&tmp,U1);
119: KSP_MatMultTranspose(ksp,Amat,U1,V1);
120: tmp = -beta; VecAXPY(&tmp,V,V1);
121: VecNorm(V1,NORM_2,&alpha);
122: tmp = 1.0 / alpha; VecScale(&tmp,V1);
124: rho = PetscSqrtScalar(rhobar*rhobar + beta*beta);
125: c = rhobar / rho;
126: s = beta / rho;
127: theta = s * alpha;
128: rhobar = - c * alpha;
129: phi = c * phibar;
130: phibar = s * phibar;
132: tmp = phi/rho;
133: VecAXPY(&tmp,W,X); /* x <- x + (phi/rho) w */
134: tmp = -theta/rho;
135: VecAYPX(&tmp,V1,W); /* w <- v - (theta/rho) w */
137: #if defined(PETSC_USE_COMPLEX)
138: rnorm = PetscRealPart(phibar);
139: #else
140: rnorm = phibar;
141: #endif
143: PetscObjectTakeAccess(ksp);
144: ksp->its++;
145: ksp->rnorm = rnorm;
146: PetscObjectGrantAccess(ksp);
147: KSPLogResidualHistory(ksp,rnorm);
148: KSPMonitor(ksp,i+1,rnorm);
149: (*ksp->converged)(ksp,i+1,rnorm,&ksp->reason,ksp->cnvP);
150: if (ksp->reason) break;
151: SWAP(U1,U,TMP);
152: SWAP(V1,V,TMP);
153: }
154: if (i == maxit) {
155: ksp->reason = KSP_DIVERGED_ITS;
156: }
158: /* KSPUnwindPreconditioner(ksp,X,W); */
160: *its = ksp->its;
161: return(0);
162: }
166: int KSPDestroy_LSQR(KSP ksp)
167: {
168: KSP_LSQR *lsqr = (KSP_LSQR*)ksp->data;
169: int ierr;
173: /* Free work vectors */
174: if (lsqr->vwork_n) {
175: VecDestroyVecs(lsqr->vwork_n,lsqr->nwork_n);
176: }
177: if (lsqr->vwork_m) {
178: VecDestroyVecs(lsqr->vwork_m,lsqr->nwork_m);
179: }
180: PetscFree(lsqr);
181: return(0);
182: }
184: EXTERN_C_BEGIN
187: int KSPCreate_LSQR(KSP ksp)
188: {
189: KSP_LSQR *lsqr;
190: int ierr;
193: PetscMalloc(sizeof(KSP_LSQR),&lsqr);
194: PetscMemzero(lsqr,sizeof(KSP_LSQR));
195: PetscLogObjectMemory(ksp,sizeof(KSP_LSQR));
196: ksp->data = (void*)lsqr;
197: ksp->pc_side = PC_LEFT;
198: ksp->ops->setup = KSPSetUp_LSQR;
199: ksp->ops->solve = KSPSolve_LSQR;
200: ksp->ops->destroy = KSPDestroy_LSQR;
201: ksp->ops->buildsolution = KSPDefaultBuildSolution;
202: ksp->ops->buildresidual = KSPDefaultBuildResidual;
203: ksp->ops->setfromoptions = 0;
204: ksp->ops->view = 0;
205: return(0);
206: }
207: EXTERN_C_END