Actual source code: symmlq.c
2: #include src/ksp/ksp/kspimpl.h
4: typedef struct {
5: PetscReal haptol;
6: } KSP_SYMMLQ;
10: PetscErrorCode KSPSetUp_SYMMLQ(KSP ksp)
11: {
15: if (ksp->pc_side == PC_RIGHT) {
16: SETERRQ(2,"No right preconditioning for KSPSYMMLQ");
17: } else if (ksp->pc_side == PC_SYMMETRIC) {
18: SETERRQ(2,"No symmetric preconditioning for KSPSYMMLQ");
19: }
20: KSPDefaultGetWork(ksp,9);
21: return(0);
22: }
26: PetscErrorCode KSPSolve_SYMMLQ(KSP ksp)
27: {
29: PetscInt 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->pc,&diagonalscale);
42: if (diagonalscale) SETERRQ1(PETSC_ERR_SUP,"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->pc,&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,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) {
79: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
80: return(0);
81: }
82: #endif
83: dp = PetscSqrtScalar(dp);
84: beta = dp; /* beta <- sqrt(r'*z) */
85: beta1 = beta;
86: s_prod = PetscAbsScalar(beta1);
88: VecCopy(R,V); /* v <- r; */
89: VecCopy(Z,U); /* u <- z; */
90: ibeta = 1.0 / beta;
91: VecScale(&ibeta,V); /* v <- ibeta*v; */
92: VecScale(&ibeta,U); /* u <- ibeta*u; */
93: VecCopy(U,Wbar); /* w_bar <- u; */
94: VecNorm(Z,NORM_2,&np); /* np <- ||z|| */
95: KSPLogResidualHistory(ksp,np);
96: KSPMonitor(ksp,0,np); /* call any registered monitor routines */
97: ksp->rnorm = np;
98: (*ksp->converged)(ksp,0,np,&ksp->reason,ksp->cnvP); /* test for convergence */
99: if (ksp->reason) return(0);
101: i = 0;
102: do {
103: ksp->its = i+1;
105: /* Update */
106: if (ksp->its > 1){
107: VecCopy(V,VOLD); /* v_old <- v; */
108: VecCopy(U,UOLD); /* u_old <- u; */
109:
110: ibeta = 1.0 / beta;
111: VecCopy(R,V);
112: VecScale(&ibeta,V); /* v <- ibeta*r; */
113: VecCopy(Z,U);
114: VecScale(&ibeta,U); /* u <- ibeta*z; */
116: VecCopy(Wbar,W);
117: VecScale(&c,W);
118: VecAXPY(&s,U,W); /* w <- c*w_bar + s*u; (w_k) */
119: ms = -s;
120: VecScale(&ms,Wbar);
121: VecAXPY(&c,U,Wbar); /* w_bar <- -s*w_bar + c*u; (w_bar_(k+1)) */
122: VecAXPY(&ceta,W,X); /* x <- x + ceta * w; (xL_k) */
124: ceta_oold = ceta_old;
125: ceta_old = ceta;
126: }
128: /* Lanczos */
129: KSP_MatMult(ksp,Amat,U,R); /* r <- Amat*u; */
130: VecDot(U,R,&alpha); /* alpha <- u'*r; */
131: KSP_PCApply(ksp,R,Z); /* z <- B*r; */
133: malpha = - alpha;
134: VecAXPY(&malpha,V,R); /* r <- r - alpha* v; */
135: VecAXPY(&malpha,U,Z); /* z <- z - alpha* u; */
136: mbeta = - beta;
137: VecAXPY(&mbeta,VOLD,R); /* r <- r - beta * v_old; */
138: VecAXPY(&mbeta,UOLD,Z); /* z <- z - beta * u_old; */
139: betaold = beta; /* beta_k */
140: VecDot(R,Z,&dp); /* dp <- r'*z; */
141: if (PetscAbsScalar(dp) < symmlq->haptol) {
142: PetscLogInfo(ksp,"KSPSolve_SYMMLQ:Detected happy breakdown %g tolerance %g\n",PetscAbsScalar(dp),symmlq->haptol);
143: dp = 0.0;
144: }
146: #if !defined(PETSC_USE_COMPLEX)
147: if (dp < 0.0) {
148: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
149: break;
150: }
151: #endif
152: beta = PetscSqrtScalar(dp); /* beta = sqrt(dp); */
154: /* QR factorization */
155: coold = cold; cold = c; soold = sold; sold = s;
156: rho0 = cold * alpha - coold * sold * betaold; /* gamma_bar */
157: rho1 = PetscSqrtScalar(rho0*rho0 + beta*beta); /* gamma */
158: rho2 = sold * alpha + coold * cold * betaold; /* delta */
159: rho3 = soold * betaold; /* epsilon */
161: /* Givens rotation: [c -s; s c] (different from the Reference!) */
162: c = rho0 / rho1; s = beta / rho1;
164: if (ksp->its==1){
165: ceta = beta1/rho1;
166: } else {
167: ceta = -(rho2*ceta_old + rho3*ceta_oold)/rho1;
168: }
169:
170: s_prod = s_prod*PetscAbsScalar(s);
171: if (c == 0.0){
172: np = s_prod*1.e16;
173: } else {
174: np = s_prod/PetscAbsScalar(c); /* residual norm for xc_k (CGNORM) */
175: }
176: ksp->rnorm = np;
177: KSPLogResidualHistory(ksp,np);
178: KSPMonitor(ksp,i+1,np);
179: (*ksp->converged)(ksp,i+1,np,&ksp->reason,ksp->cnvP); /* test for convergence */
180: if (ksp->reason) break;
181: i++;
182: } while (i<ksp->max_it);
184: /* move to the CG point: xc_(k+1) */
185: if (c == 0.0){
186: ceta_bar = ceta*1.e15;
187: } else {
188: ceta_bar = ceta/c;
189: }
190: VecAXPY(&ceta_bar,Wbar,X); /* x <- x + ceta_bar*w_bar */
192: if (i >= ksp->max_it) {
193: ksp->reason = KSP_DIVERGED_ITS;
194: }
195: return(0);
196: }
198: /*MC
199: KSPSYMMLQ - This code implements the SYMMLQ method.
201: Options Database Keys:
202: . see KSPSolve()
204: Level: beginner
206: Notes: Reference: Paige & Saunders, 1975.
208: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP
209: M*/
213: PetscErrorCode KSPCreate_SYMMLQ(KSP ksp)
214: {
215: KSP_SYMMLQ *symmlq;
219: ksp->pc_side = PC_LEFT;
221: PetscNew(KSP_SYMMLQ,&symmlq);
222: symmlq->haptol = 1.e-18;
223: ksp->data = (void*)symmlq;
225: /*
226: Sets the functions that are associated with this data structure
227: (in C++ this is the same as defining virtual functions)
228: */
229: ksp->ops->setup = KSPSetUp_SYMMLQ;
230: ksp->ops->solve = KSPSolve_SYMMLQ;
231: ksp->ops->destroy = KSPDefaultDestroy;
232: ksp->ops->setfromoptions = 0;
233: ksp->ops->buildsolution = KSPDefaultBuildSolution;
234: ksp->ops->buildresidual = KSPDefaultBuildResidual;
235: return(0);
236: }