Actual source code: symmlq.c
1: /*$Id: symmlq.c,v 1.13 2001/03/23 23:23:49 balay Exp $*/
2: /*
3: This code implements the SYMMLQ method.
4: Reference: Paige & Saunders, 1975.
5: */
6: #include "src/sles/ksp/kspimpl.h"
8: typedef struct {
9: double haptol;
10: } KSP_SYMMLQ;
12: int KSPSetUp_SYMMLQ(KSP ksp)
13: {
17: if (ksp->pc_side == PC_RIGHT) {
18: SETERRQ(2,"No right preconditioning for KSPSYMMLQ");
19: } else if (ksp->pc_side == PC_SYMMETRIC) {
20: SETERRQ(2,"No symmetric preconditioning for KSPSYMMLQ");
21: }
22: KSPDefaultGetWork(ksp,9);
23: return(0);
24: }
26: int KSPSolve_SYMMLQ(KSP ksp,int *its)
27: {
28: int ierr,i,maxit;
29: Scalar alpha,malpha,beta,mbeta,ibeta,betaold,beta1,ceta,ceta_oold = 0.0, ceta_old = 0.0,ceta_bar;
30: Scalar c=1.0,cold=1.0,s=0.0,sold=0.0,coold,soold,ms,rho0,rho1,rho2,rho3;
31: Scalar mone = -1.0,zero = 0.0,dp = 0.0;
32: PetscReal np,s_prod;
33: Vec X,B,R,Z,U,V,W,UOLD,VOLD,Wbar;
34: Mat Amat,Pmat;
35: MatStructure pflag;
36: KSP_SYMMLQ *symmlq = (KSP_SYMMLQ*)ksp->data;
37: PetscTruth diagonalscale;
40: ierr = PCDiagonalScale(ksp->B,&diagonalscale);
41: if (diagonalscale) SETERRQ1(1,"Krylov method %s does not support diagonal scaling",ksp->type_name);
43: maxit = ksp->max_it;
44: X = ksp->vec_sol;
45: B = ksp->vec_rhs;
46: R = ksp->work[0];
47: Z = ksp->work[1];
48: U = ksp->work[2];
49: V = ksp->work[3];
50: W = ksp->work[4];
51: UOLD = ksp->work[5];
52: VOLD = ksp->work[6];
53: Wbar = ksp->work[7];
54:
55: PCGetOperators(ksp->B,&Amat,&Pmat,&pflag);
57: ksp->its = 0;
59: VecSet(&zero,UOLD); /* u_old <- zeros; */
60: VecCopy(UOLD,VOLD); /* v_old <- u_old; */
61: VecCopy(UOLD,W); /* w <- u_old; */
62: VecCopy(UOLD,Wbar); /* w_bar <- u_old; */
63: if (!ksp->guess_zero) {
64: KSP_MatMult(ksp,Amat,X,R); /* r <- b - A*x */
65: VecAYPX(&mone,B,R);
66: } else {
67: VecCopy(B,R); /* r <- b (x is 0) */
68: }
70: KSP_PCApply(ksp,ksp->B,R,Z); /* z <- B*r */
71: VecDot(R,Z,&dp); /* dp = r'*z; */
72: if (PetscAbsScalar(dp) < symmlq->haptol) {
73: PetscLogInfo(ksp,"KSPSolve_SYMMLQ:Detected happy breakdown %g tolerance %gn",PetscAbsScalar(dp),symmlq->haptol);
74: dp = 0.0;
75: }
77: #if !defined(PETSC_USE_COMPLEX)
78: if (dp < 0.0) SETERRQ(PETSC_ERR_KSP_BRKDWN,"Indefinite preconditioner");
79: #endif
80: dp = PetscSqrtScalar(dp);
81: beta = dp; /* beta <- sqrt(r'*z) */
82: beta1 = beta;
83: s_prod = PetscAbsScalar(beta1);
85: VecCopy(R,V); /* v <- r; */
86: VecCopy(Z,U); /* u <- z; */
87: ibeta = 1.0 / beta;
88: VecScale(&ibeta,V); /* v <- ibeta*v; */
89: VecScale(&ibeta,U); /* u <- ibeta*u; */
90: VecCopy(U,Wbar); /* w_bar <- u; */
91: VecNorm(Z,NORM_2,&np); /* np <- ||z|| */
92: (*ksp->converged)(ksp,0,np,&ksp->reason,ksp->cnvP); /* test for convergence */
93: if (ksp->reason) {*its = 0; return(0);}
94: KSPLogResidualHistory(ksp,np);
95: KSPMonitor(ksp,0,np); /* call any registered monitor routines */
96: ksp->rnorm = np;
98: for (i=0; i<maxit; i++){
99: ksp->its = i+1;
101: /* Update */
102: if (ksp->its > 1){
103: VecCopy(V,VOLD); /* v_old <- v; */
104: VecCopy(U,UOLD); /* u_old <- u; */
105:
106: ibeta = 1.0 / beta;
107: VecCopy(R,V);
108: VecScale(&ibeta,V); /* v <- ibeta*r; */
109: VecCopy(Z,U);
110: VecScale(&ibeta,U); /* u <- ibeta*z; */
112: VecCopy(Wbar,W);
113: VecScale(&c,W);
114: VecAXPY(&s,U,W); /* w <- c*w_bar + s*u; (w_k) */
115: ms = -s;
116: VecScale(&ms,Wbar);
117: VecAXPY(&c,U,Wbar); /* w_bar <- -s*w_bar + c*u; (w_bar_(k+1)) */
118: VecAXPY(&ceta,W,X); /* x <- x + ceta * w; (xL_k) */
120: ceta_oold = ceta_old;
121: ceta_old = ceta;
122: }
124: /* Lanczos */
125: KSP_MatMult(ksp,Amat,U,R); /* r <- Amat*u; */
126: VecDot(U,R,&alpha); /* alpha <- u'*r; */
127: KSP_PCApply(ksp,ksp->B,R,Z); /* z <- B*r; */
129: malpha = - alpha;
130: VecAXPY(&malpha,V,R); /* r <- r - alpha* v; */
131: VecAXPY(&malpha,U,Z); /* z <- z - alpha* u; */
132: mbeta = - beta;
133: VecAXPY(&mbeta,VOLD,R); /* r <- r - beta * v_old; */
134: VecAXPY(&mbeta,UOLD,Z); /* z <- z - beta * u_old; */
135: betaold = beta; /* beta_k */
136: VecDot(R,Z,&dp); /* dp <- r'*z; */
137: if (PetscAbsScalar(dp) < symmlq->haptol) {
138: PetscLogInfo(ksp,"KSPSolve_SYMMLQ:Detected happy breakdown %g tolerance %gn",PetscAbsScalar(dp),symmlq->haptol);
139: dp = 0.0;
140: }
142: #if !defined(PETSC_USE_COMPLEX)
143: if (dp < 0.0) SETERRQ(PETSC_ERR_KSP_BRKDWN,"Indefinite preconditioner");
144: #endif
145: beta = PetscSqrtScalar(dp); /* beta = sqrt(dp); */
147: /* QR factorization */
148: coold = cold; cold = c; soold = sold; sold = s;
149: rho0 = cold * alpha - coold * sold * betaold; /* gamma_bar */
150: rho1 = PetscSqrtScalar(rho0*rho0 + beta*beta); /* gamma */
151: rho2 = sold * alpha + coold * cold * betaold; /* delta */
152: rho3 = soold * betaold; /* epsilon */
154: /* Givens rotation: [c -s; s c] (different from the Reference!) */
155: c = rho0 / rho1; s = beta / rho1;
157: if (ksp->its==1){
158: ceta = beta1/rho1;
159: } else {
160: ceta = -(rho2*ceta_old + rho3*ceta_oold)/rho1;
161: }
162:
163: s_prod = s_prod*PetscAbsScalar(s);
164: if (c == 0.0){
165: np = s_prod*1.e16;
166: } else {
167: np = s_prod/PetscAbsScalar(c); /* residual norm for xc_k (CGNORM) */
168: }
169: ksp->rnorm = np;
170: KSPLogResidualHistory(ksp,np);
171: KSPMonitor(ksp,i+1,np);
172: (*ksp->converged)(ksp,i+1,np,&ksp->reason,ksp->cnvP); /* test for convergence */
173: if (ksp->reason) break;
174: }
176: /* move to the CG point: xc_(k+1) */
177: if (c == 0.0){
178: ceta_bar = ceta*1.e15;
179: } else {
180: ceta_bar = ceta/c;
181: }
182: VecAXPY(&ceta_bar,Wbar,X); /* x <- x + ceta_bar*w_bar */
184: if (i == maxit) {
185: ksp->reason = KSP_DIVERGED_ITS;
186: }
187: *its = ksp->its;
189: return(0);
190: }
192: EXTERN_C_BEGIN
193: int KSPCreate_SYMMLQ(KSP ksp)
194: {
195: KSP_SYMMLQ *symmlq;
200: ksp->pc_side = PC_LEFT;
201: ksp->calc_res = PETSC_TRUE;
203: ierr = PetscNew(KSP_SYMMLQ,&symmlq);
204: symmlq->haptol = 1.e-18;
205: ksp->data = (void*)symmlq;
207: /*
208: Sets the functions that are associated with this data structure
209: (in C++ this is the same as defining virtual functions)
210: */
211: ksp->ops->setup = KSPSetUp_SYMMLQ;
212: ksp->ops->solve = KSPSolve_SYMMLQ;
213: ksp->ops->destroy = KSPDefaultDestroy;
214: ksp->ops->setfromoptions = 0;
215: ksp->ops->buildsolution = KSPDefaultBuildSolution;
216: ksp->ops->buildresidual = KSPDefaultBuildResidual;
218: return(0);
219: }
220: EXTERN_C_END