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: #include src/ksp/ksp/kspimpl.h
7: typedef struct {
8: int 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 m */
11: } KSP_LSQR;
15: static int KSPSetUp_LSQR(KSP ksp)
16: {
17: int ierr,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: VecDuplicateVecs(ksp->vec_rhs,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: VecDuplicateVecs(ksp->vec_sol,nw,&lsqr->vwork_n);
38: PetscLogObjectParents(ksp,nw,lsqr->vwork_n);
40: return(0);
41: }
45: static int KSPSolve_LSQR(KSP ksp)
46: {
47: int i,ierr;
48: PetscScalar rho,rhobar,phi,phibar,theta,c,s,tmp,zero = 0.0,mone=-1.0;
49: PetscReal beta,alpha,rnorm;
50: Vec X,B,V,V1,U,U1,TMP,W;
51: Mat Amat,Pmat;
52: MatStructure pflag;
53: KSP_LSQR *lsqr = (KSP_LSQR*)ksp->data;
54: PetscTruth diagonalscale;
57: PCDiagonalScale(ksp->B,&diagonalscale);
58: if (diagonalscale) SETERRQ1(1,"Krylov method %s does not support diagonal scaling",ksp->type_name);
60: PCGetOperators(ksp->B,&Amat,&Pmat,&pflag);
62: /* vectors of length m, where system size is mxn */
63: B = ksp->vec_rhs;
64: U = lsqr->vwork_m[0];
65: U1 = lsqr->vwork_m[1];
67: /* vectors of length n */
68: X = ksp->vec_sol;
69: W = lsqr->vwork_n[0];
70: V = lsqr->vwork_n[1];
71: V1 = lsqr->vwork_n[2];
73: /* Compute initial residual, temporarily use work vector u */
74: if (!ksp->guess_zero) {
75: KSP_MatMult(ksp,Amat,X,U); /* u <- b - Ax */
76: VecAYPX(&mone,B,U);
77: } else {
78: VecCopy(B,U); /* u <- b (x is 0) */
79: }
81: /* Test for nothing to do */
82: VecNorm(U,NORM_2,&rnorm);
83: PetscObjectTakeAccess(ksp);
84: ksp->its = 0;
85: ksp->rnorm = rnorm;
86: PetscObjectGrantAccess(ksp);
87: KSPLogResidualHistory(ksp,rnorm);
88: KSPMonitor(ksp,0,rnorm);
89: (*ksp->converged)(ksp,0,rnorm,&ksp->reason,ksp->cnvP);
90: if (ksp->reason) return(0);
92: VecCopy(B,U);
93: VecNorm(U,NORM_2,&beta);
94: tmp = 1.0/beta; VecScale(&tmp,U);
95: KSP_MatMultTranspose(ksp,Amat,U,V);
96: VecNorm(V,NORM_2,&alpha);
97: tmp = 1.0/alpha; VecScale(&tmp,V);
99: VecCopy(V,W);
100: VecSet(&zero,X);
102: phibar = beta;
103: rhobar = alpha;
104: i = 0;
105: do {
107: KSP_MatMult(ksp,Amat,V,U1);
108: tmp = -alpha; VecAXPY(&tmp,U,U1);
109: VecNorm(U1,NORM_2,&beta);
110: tmp = 1.0/beta; VecScale(&tmp,U1);
112: KSP_MatMultTranspose(ksp,Amat,U1,V1);
113: tmp = -beta; VecAXPY(&tmp,V,V1);
114: VecNorm(V1,NORM_2,&alpha);
115: tmp = 1.0 / alpha; VecScale(&tmp,V1);
117: rho = PetscSqrtScalar(rhobar*rhobar + beta*beta);
118: c = rhobar / rho;
119: s = beta / rho;
120: theta = s * alpha;
121: rhobar = - c * alpha;
122: phi = c * phibar;
123: phibar = s * phibar;
125: tmp = phi/rho;
126: VecAXPY(&tmp,W,X); /* x <- x + (phi/rho) w */
127: tmp = -theta/rho;
128: VecAYPX(&tmp,V1,W); /* w <- v - (theta/rho) w */
130: #if defined(PETSC_USE_COMPLEX)
131: rnorm = PetscRealPart(phibar);
132: #else
133: rnorm = phibar;
134: #endif
136: PetscObjectTakeAccess(ksp);
137: ksp->its++;
138: ksp->rnorm = rnorm;
139: PetscObjectGrantAccess(ksp);
140: KSPLogResidualHistory(ksp,rnorm);
141: KSPMonitor(ksp,i+1,rnorm);
142: (*ksp->converged)(ksp,i+1,rnorm,&ksp->reason,ksp->cnvP);
143: if (ksp->reason) break;
144: SWAP(U1,U,TMP);
145: SWAP(V1,V,TMP);
147: i++;
148: } while (i<ksp->max_it);
149: if (i == ksp->max_it) {
150: ksp->reason = KSP_DIVERGED_ITS;
151: }
153: /* KSPUnwindPreconditioner(ksp,X,W); */
155: return(0);
156: }
160: int KSPDestroy_LSQR(KSP ksp)
161: {
162: KSP_LSQR *lsqr = (KSP_LSQR*)ksp->data;
163: int ierr;
167: /* Free work vectors */
168: if (lsqr->vwork_n) {
169: VecDestroyVecs(lsqr->vwork_n,lsqr->nwork_n);
170: }
171: if (lsqr->vwork_m) {
172: VecDestroyVecs(lsqr->vwork_m,lsqr->nwork_m);
173: }
174: PetscFree(lsqr);
175: return(0);
176: }
178: /*MC
179: KSPLSQR - This implements LSQR (Paige and Saunders, ACM Transactions on
180: Mathematical Software, Vol 8, pp 43-71, 1982).
182: Options Database Keys:
183: . see KSPSolve()
185: Level: beginner
187: Notes: This algorithm DOES NOT use a preconditioner. It ignores
188: any preconditioner arguments specified.
190: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP
192: M*/
193: EXTERN_C_BEGIN
196: int KSPCreate_LSQR(KSP ksp)
197: {
198: KSP_LSQR *lsqr;
199: int ierr;
202: PetscMalloc(sizeof(KSP_LSQR),&lsqr);
203: PetscMemzero(lsqr,sizeof(KSP_LSQR));
204: PetscLogObjectMemory(ksp,sizeof(KSP_LSQR));
205: ksp->data = (void*)lsqr;
206: ksp->pc_side = PC_LEFT;
207: ksp->ops->setup = KSPSetUp_LSQR;
208: ksp->ops->solve = KSPSolve_LSQR;
209: ksp->ops->destroy = KSPDestroy_LSQR;
210: ksp->ops->buildsolution = KSPDefaultBuildSolution;
211: ksp->ops->buildresidual = KSPDefaultBuildResidual;
212: ksp->ops->setfromoptions = 0;
213: ksp->ops->view = 0;
214: return(0);
215: }
216: EXTERN_C_END