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;
20: static int KSPSetUp_LSQR(KSP ksp)
21: {
22: int ierr,nw;
23: KSP_LSQR *lsqr = (KSP_LSQR*)ksp->data;
26: if (ksp->pc_side == PC_SYMMETRIC){
27: SETERRQ(2,"no symmetric preconditioning for KSPLSQR");
28: }
30: /* Get work vectors */
31: lsqr->nwork_m = nw = 2;
32: if (lsqr->vwork_m) {
33: VecDestroyVecs(lsqr->vwork_m,lsqr->nwork_m);
34: }
35: VecDuplicateVecs(ksp->vec_rhs,nw,&lsqr->vwork_m);
36: PetscLogObjectParents(ksp,nw,lsqr->vwork_m);
38: lsqr->nwork_n = nw = 3;
39: if (lsqr->vwork_n) {
40: VecDestroyVecs(lsqr->vwork_n,lsqr->nwork_n);
41: }
42: VecDuplicateVecs(ksp->vec_sol,nw,&lsqr->vwork_n);
43: PetscLogObjectParents(ksp,nw,lsqr->vwork_n);
45: return(0);
46: }
48: static int KSPSolve_LSQR(KSP ksp,int *its)
49: {
50: int i,maxit,ierr;
51: PetscScalar rho,rhobar,phi,phibar,theta,c,s,tmp,zero = 0.0,mone=-1.0;
52: PetscReal beta,alpha,rnorm;
53: Vec X,B,V,V1,U,U1,TMP,W;
54: Mat Amat,Pmat;
55: MatStructure pflag;
56: KSP_LSQR *lsqr = (KSP_LSQR*)ksp->data;
57: PetscTruth diagonalscale;
60: ierr = PCDiagonalScale(ksp->B,&diagonalscale);
61: if (diagonalscale) SETERRQ1(1,"Krylov method %s does not support diagonal scaling",ksp->type_name);
63: ierr = PCGetOperators(ksp->B,&Amat,&Pmat,&pflag);
64: maxit = ksp->max_it;
66: /* vectors of length m, where system size is mxn */
67: B = ksp->vec_rhs;
68: U = lsqr->vwork_m[0];
69: U1 = lsqr->vwork_m[1];
71: /* vectors of length n */
72: X = ksp->vec_sol;
73: W = lsqr->vwork_n[0];
74: V = lsqr->vwork_n[1];
75: V1 = lsqr->vwork_n[2];
77: /* Compute initial residual, temporarily use work vector u */
78: if (!ksp->guess_zero) {
79: KSP_MatMult(ksp,Amat,X,U); /* u <- b - Ax */
80: VecAYPX(&mone,B,U);
81: } else {
82: VecCopy(B,U); /* u <- b (x is 0) */
83: }
85: /* Test for nothing to do */
86: VecNorm(U,NORM_2,&rnorm);
87: PetscObjectTakeAccess(ksp);
88: ksp->its = 0;
89: ksp->rnorm = rnorm;
90: PetscObjectGrantAccess(ksp);
91: (*ksp->converged)(ksp,0,rnorm,&ksp->reason,ksp->cnvP);
92: if (ksp->reason) {*its = 0; return(0);}
93: KSPLogResidualHistory(ksp,rnorm);
94: KSPMonitor(ksp,0,rnorm);
96: VecCopy(B,U);
97: VecNorm(U,NORM_2,&beta);
98: tmp = 1.0/beta; VecScale(&tmp,U);
99: KSP_MatMultTranspose(ksp,Amat,U,V);
100: VecNorm(V,NORM_2,&alpha);
101: tmp = 1.0/alpha; VecScale(&tmp,V);
103: VecCopy(V,W);
104: VecSet(&zero,X);
106: phibar = beta;
107: rhobar = alpha;
108: for (i=0; i<maxit; i++) {
110: KSP_MatMult(ksp,Amat,V,U1);
111: tmp = -alpha; VecAXPY(&tmp,U,U1);
112: VecNorm(U1,NORM_2,&beta);
113: tmp = 1.0/beta; VecScale(&tmp,U1);
115: KSP_MatMultTranspose(ksp,Amat,U1,V1);
116: tmp = -beta; VecAXPY(&tmp,V,V1);
117: VecNorm(V1,NORM_2,&alpha);
118: tmp = 1.0 / alpha; VecScale(&tmp,V1);
120: rho = PetscSqrtScalar(rhobar*rhobar + beta*beta);
121: c = rhobar / rho;
122: s = beta / rho;
123: theta = s * alpha;
124: rhobar = - c * alpha;
125: phi = c * phibar;
126: phibar = s * phibar;
128: tmp = phi/rho;
129: VecAXPY(&tmp,W,X); /* x <- x + (phi/rho) w */
130: tmp = -theta/rho;
131: VecAYPX(&tmp,V1,W); /* w <- v - (theta/rho) w */
133: #if defined(PETSC_USE_COMPLEX)
134: rnorm = PetscRealPart(phibar);
135: #else
136: rnorm = phibar;
137: #endif
139: PetscObjectTakeAccess(ksp);
140: ksp->its++;
141: ksp->rnorm = rnorm;
142: PetscObjectGrantAccess(ksp);
143: KSPLogResidualHistory(ksp,rnorm);
144: KSPMonitor(ksp,i+1,rnorm);
145: (*ksp->converged)(ksp,i+1,rnorm,&ksp->reason,ksp->cnvP);
146: if (ksp->reason) break;
147: SWAP(U1,U,TMP);
148: SWAP(V1,V,TMP);
149: }
150: if (i == maxit) {
151: i--;
152: ksp->reason = KSP_DIVERGED_ITS;
153: }
155: /* KSPUnwindPreconditioner(ksp,X,W); */
157: *its = i + 1;
158: return(0);
159: }
161: int KSPDestroy_LSQR(KSP ksp)
162: {
163: KSP_LSQR *lsqr = (KSP_LSQR*)ksp->data;
164: int ierr;
168: /* Free work vectors */
169: if (lsqr->vwork_n) {
170: VecDestroyVecs(lsqr->vwork_n,lsqr->nwork_n);
171: }
172: if (lsqr->vwork_m) {
173: VecDestroyVecs(lsqr->vwork_m,lsqr->nwork_m);
174: }
175: PetscFree(lsqr);
176: return(0);
177: }
179: EXTERN_C_BEGIN
180: int KSPCreate_LSQR(KSP ksp)
181: {
182: KSP_LSQR *lsqr;
183: int ierr;
186: PetscMalloc(sizeof(KSP_LSQR),&lsqr);
187: PetscMemzero(lsqr,sizeof(KSP_LSQR));
188: PetscLogObjectMemory(ksp,sizeof(KSP_LSQR));
189: ksp->data = (void*)lsqr;
190: ksp->pc_side = PC_LEFT;
191: ksp->ops->setup = KSPSetUp_LSQR;
192: ksp->ops->solve = KSPSolve_LSQR;
193: ksp->ops->destroy = KSPDestroy_LSQR;
194: ksp->ops->buildsolution = KSPDefaultBuildSolution;
195: ksp->ops->buildresidual = KSPDefaultBuildResidual;
196: ksp->ops->setfromoptions = 0;
197: ksp->ops->view = 0;
198: return(0);
199: }
200: EXTERN_C_END