Actual source code: symmlq.c
1: /*$Id: symmlq.c,v 1.16 2001/08/07 03:03:56 balay Exp $*/
3: #include src/ksp/ksp/kspimpl.h
5: typedef struct {
6: PetscReal haptol;
7: } KSP_SYMMLQ;
11: int KSPSetUp_SYMMLQ(KSP ksp)
12: {
16: if (ksp->pc_side == PC_RIGHT) {
17: SETERRQ(2,"No right preconditioning for KSPSYMMLQ");
18: } else if (ksp->pc_side == PC_SYMMETRIC) {
19: SETERRQ(2,"No symmetric preconditioning for KSPSYMMLQ");
20: }
21: KSPDefaultGetWork(ksp,9);
22: return(0);
23: }
27: int KSPSolve_SYMMLQ(KSP ksp)
28: {
29: int ierr,i;
30: PetscScalar alpha,malpha,beta,mbeta,ibeta,betaold,beta1,ceta,ceta_oold = 0.0, ceta_old = 0.0,ceta_bar;
31: PetscScalar c=1.0,cold=1.0,s=0.0,sold=0.0,coold,soold,ms,rho0,rho1,rho2,rho3;
32: PetscScalar mone = -1.0,zero = 0.0,dp = 0.0;
33: PetscReal np,s_prod;
34: Vec X,B,R,Z,U,V,W,UOLD,VOLD,Wbar;
35: Mat Amat,Pmat;
36: MatStructure pflag;
37: KSP_SYMMLQ *symmlq = (KSP_SYMMLQ*)ksp->data;
38: PetscTruth diagonalscale;
41: PCDiagonalScale(ksp->B,&diagonalscale);
42: if (diagonalscale) SETERRQ1(1,"Krylov method %s does not support diagonal scaling",ksp->type_name);
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 %g\n",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: KSPLogResidualHistory(ksp,np);
93: KSPMonitor(ksp,0,np); /* call any registered monitor routines */
94: ksp->rnorm = np;
95: (*ksp->converged)(ksp,0,np,&ksp->reason,ksp->cnvP); /* test for convergence */
96: if (ksp->reason) return(0);
98: i = 0;
99: do {
100: ksp->its = i+1;
102: /* Update */
103: if (ksp->its > 1){
104: VecCopy(V,VOLD); /* v_old <- v; */
105: VecCopy(U,UOLD); /* u_old <- u; */
106:
107: ibeta = 1.0 / beta;
108: VecCopy(R,V);
109: VecScale(&ibeta,V); /* v <- ibeta*r; */
110: VecCopy(Z,U);
111: VecScale(&ibeta,U); /* u <- ibeta*z; */
113: VecCopy(Wbar,W);
114: VecScale(&c,W);
115: VecAXPY(&s,U,W); /* w <- c*w_bar + s*u; (w_k) */
116: ms = -s;
117: VecScale(&ms,Wbar);
118: VecAXPY(&c,U,Wbar); /* w_bar <- -s*w_bar + c*u; (w_bar_(k+1)) */
119: VecAXPY(&ceta,W,X); /* x <- x + ceta * w; (xL_k) */
121: ceta_oold = ceta_old;
122: ceta_old = ceta;
123: }
125: /* Lanczos */
126: KSP_MatMult(ksp,Amat,U,R); /* r <- Amat*u; */
127: VecDot(U,R,&alpha); /* alpha <- u'*r; */
128: KSP_PCApply(ksp,ksp->B,R,Z); /* z <- B*r; */
130: malpha = - alpha;
131: VecAXPY(&malpha,V,R); /* r <- r - alpha* v; */
132: VecAXPY(&malpha,U,Z); /* z <- z - alpha* u; */
133: mbeta = - beta;
134: VecAXPY(&mbeta,VOLD,R); /* r <- r - beta * v_old; */
135: VecAXPY(&mbeta,UOLD,Z); /* z <- z - beta * u_old; */
136: betaold = beta; /* beta_k */
137: VecDot(R,Z,&dp); /* dp <- r'*z; */
138: if (PetscAbsScalar(dp) < symmlq->haptol) {
139: PetscLogInfo(ksp,"KSPSolve_SYMMLQ:Detected happy breakdown %g tolerance %g\n",PetscAbsScalar(dp),symmlq->haptol);
140: dp = 0.0;
141: }
143: #if !defined(PETSC_USE_COMPLEX)
144: if (dp < 0.0) SETERRQ(PETSC_ERR_KSP_BRKDWN,"Indefinite preconditioner");
145: #endif
146: beta = PetscSqrtScalar(dp); /* beta = sqrt(dp); */
148: /* QR factorization */
149: coold = cold; cold = c; soold = sold; sold = s;
150: rho0 = cold * alpha - coold * sold * betaold; /* gamma_bar */
151: rho1 = PetscSqrtScalar(rho0*rho0 + beta*beta); /* gamma */
152: rho2 = sold * alpha + coold * cold * betaold; /* delta */
153: rho3 = soold * betaold; /* epsilon */
155: /* Givens rotation: [c -s; s c] (different from the Reference!) */
156: c = rho0 / rho1; s = beta / rho1;
158: if (ksp->its==1){
159: ceta = beta1/rho1;
160: } else {
161: ceta = -(rho2*ceta_old + rho3*ceta_oold)/rho1;
162: }
163:
164: s_prod = s_prod*PetscAbsScalar(s);
165: if (c == 0.0){
166: np = s_prod*1.e16;
167: } else {
168: np = s_prod/PetscAbsScalar(c); /* residual norm for xc_k (CGNORM) */
169: }
170: ksp->rnorm = np;
171: KSPLogResidualHistory(ksp,np);
172: KSPMonitor(ksp,i+1,np);
173: (*ksp->converged)(ksp,i+1,np,&ksp->reason,ksp->cnvP); /* test for convergence */
174: if (ksp->reason) break;
175: i++;
176: } while (i<ksp->max_it);
178: /* move to the CG point: xc_(k+1) */
179: if (c == 0.0){
180: ceta_bar = ceta*1.e15;
181: } else {
182: ceta_bar = ceta/c;
183: }
184: VecAXPY(&ceta_bar,Wbar,X); /* x <- x + ceta_bar*w_bar */
186: if (i == ksp->max_it) {
187: ksp->reason = KSP_DIVERGED_ITS;
188: }
189: return(0);
190: }
192: /*MC
193: KSPSYMMLQ - This code implements the SYMMLQ method.
194: Reference: Paige & Saunders, 1975.
196: Options Database Keys:
197: . see KSPSolve()
199: Level: beginner
201: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP
203: M*/
205: EXTERN_C_BEGIN
208: int KSPCreate_SYMMLQ(KSP ksp)
209: {
210: KSP_SYMMLQ *symmlq;
215: ksp->pc_side = PC_LEFT;
217: PetscNew(KSP_SYMMLQ,&symmlq);
218: symmlq->haptol = 1.e-18;
219: ksp->data = (void*)symmlq;
221: /*
222: Sets the functions that are associated with this data structure
223: (in C++ this is the same as defining virtual functions)
224: */
225: ksp->ops->setup = KSPSetUp_SYMMLQ;
226: ksp->ops->solve = KSPSolve_SYMMLQ;
227: ksp->ops->destroy = KSPDefaultDestroy;
228: ksp->ops->setfromoptions = 0;
229: ksp->ops->buildsolution = KSPDefaultBuildSolution;
230: ksp->ops->buildresidual = KSPDefaultBuildResidual;
232: return(0);
233: }
234: EXTERN_C_END