Actual source code: symmlq.c
2: #include <petsc/private/kspimpl.h>
4: typedef struct {
5: PetscReal haptol;
6: } KSP_SYMMLQ;
8: PetscErrorCode KSPSetUp_SYMMLQ(KSP ksp)
9: {
10: KSPSetWorkVecs(ksp, 9);
11: return 0;
12: }
14: PetscErrorCode KSPSolve_SYMMLQ(KSP ksp)
15: {
16: PetscInt i;
17: PetscScalar alpha, beta, ibeta, betaold, beta1, ceta = 0, ceta_oold = 0.0, ceta_old = 0.0, ceta_bar;
18: PetscScalar c = 1.0, cold = 1.0, s = 0.0, sold = 0.0, coold, soold, rho0, rho1, rho2, rho3;
19: PetscScalar dp = 0.0;
20: PetscReal np = 0.0, s_prod;
21: Vec X, B, R, Z, U, V, W, UOLD, VOLD, Wbar;
22: Mat Amat, Pmat;
23: KSP_SYMMLQ *symmlq = (KSP_SYMMLQ *)ksp->data;
24: PetscBool diagonalscale;
26: PCGetDiagonalScale(ksp->pc, &diagonalscale);
29: X = ksp->vec_sol;
30: B = ksp->vec_rhs;
31: R = ksp->work[0];
32: Z = ksp->work[1];
33: U = ksp->work[2];
34: V = ksp->work[3];
35: W = ksp->work[4];
36: UOLD = ksp->work[5];
37: VOLD = ksp->work[6];
38: Wbar = ksp->work[7];
40: PCGetOperators(ksp->pc, &Amat, &Pmat);
42: ksp->its = 0;
44: VecSet(UOLD, 0.0); /* u_old <- zeros; */
45: VecCopy(UOLD, VOLD); /* v_old <- u_old; */
46: VecCopy(UOLD, W); /* w <- u_old; */
47: VecCopy(UOLD, Wbar); /* w_bar <- u_old; */
48: if (!ksp->guess_zero) {
49: KSP_MatMult(ksp, Amat, X, R); /* r <- b - A*x */
50: VecAYPX(R, -1.0, B);
51: } else {
52: VecCopy(B, R); /* r <- b (x is 0) */
53: }
55: KSP_PCApply(ksp, R, Z); /* z <- B*r */
56: VecDot(R, Z, &dp); /* dp = r'*z; */
57: KSPCheckDot(ksp, dp);
58: if (PetscAbsScalar(dp) < symmlq->haptol) {
59: PetscInfo(ksp, "Detected happy breakdown %g tolerance %g\n", (double)PetscAbsScalar(dp), (double)symmlq->haptol);
60: ksp->rnorm = 0.0; /* what should we really put here? */
61: ksp->reason = KSP_CONVERGED_HAPPY_BREAKDOWN; /* bugfix proposed by Lourens (lourens.vanzanen@shell.com) */
62: return 0;
63: }
65: #if !defined(PETSC_USE_COMPLEX)
66: if (dp < 0.0) {
67: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
68: return 0;
69: }
70: #endif
71: dp = PetscSqrtScalar(dp);
72: beta = dp; /* beta <- sqrt(r'*z) */
73: beta1 = beta;
74: s_prod = PetscAbsScalar(beta1);
76: VecCopy(R, V); /* v <- r; */
77: VecCopy(Z, U); /* u <- z; */
78: ibeta = 1.0 / beta;
79: VecScale(V, ibeta); /* v <- ibeta*v; */
80: VecScale(U, ibeta); /* u <- ibeta*u; */
81: VecCopy(U, Wbar); /* w_bar <- u; */
82: if (ksp->normtype != KSP_NORM_NONE) {
83: VecNorm(Z, NORM_2, &np); /* np <- ||z|| */
84: KSPCheckNorm(ksp, np);
85: }
86: KSPLogResidualHistory(ksp, np);
87: KSPMonitor(ksp, 0, np);
88: ksp->rnorm = np;
89: (*ksp->converged)(ksp, 0, np, &ksp->reason, ksp->cnvP); /* test for convergence */
90: if (ksp->reason) return 0;
92: i = 0;
93: ceta = 0.;
94: do {
95: ksp->its = i + 1;
97: /* Update */
98: if (ksp->its > 1) {
99: VecCopy(V, VOLD); /* v_old <- v; */
100: VecCopy(U, UOLD); /* u_old <- u; */
102: VecCopy(R, V);
103: VecScale(V, 1.0 / beta); /* v <- ibeta*r; */
104: VecCopy(Z, U);
105: VecScale(U, 1.0 / beta); /* u <- ibeta*z; */
107: VecCopy(Wbar, W);
108: VecScale(W, c);
109: VecAXPY(W, s, U); /* w <- c*w_bar + s*u; (w_k) */
110: VecScale(Wbar, -s);
111: VecAXPY(Wbar, c, U); /* w_bar <- -s*w_bar + c*u; (w_bar_(k+1)) */
112: VecAXPY(X, ceta, W); /* x <- x + ceta * w; (xL_k) */
114: ceta_oold = ceta_old;
115: ceta_old = ceta;
116: }
118: /* Lanczos */
119: KSP_MatMult(ksp, Amat, U, R); /* r <- Amat*u; */
120: VecDot(U, R, &alpha); /* alpha <- u'*r; */
121: KSP_PCApply(ksp, R, Z); /* z <- B*r; */
123: VecAXPY(R, -alpha, V); /* r <- r - alpha* v; */
124: VecAXPY(Z, -alpha, U); /* z <- z - alpha* u; */
125: VecAXPY(R, -beta, VOLD); /* r <- r - beta * v_old; */
126: VecAXPY(Z, -beta, UOLD); /* z <- z - beta * u_old; */
127: betaold = beta; /* beta_k */
128: VecDot(R, Z, &dp); /* dp <- r'*z; */
129: KSPCheckDot(ksp, dp);
130: if (PetscAbsScalar(dp) < symmlq->haptol) {
131: PetscInfo(ksp, "Detected happy breakdown %g tolerance %g\n", (double)PetscAbsScalar(dp), (double)symmlq->haptol);
132: dp = 0.0;
133: }
135: #if !defined(PETSC_USE_COMPLEX)
136: if (dp < 0.0) {
137: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
138: break;
139: }
140: #endif
141: beta = PetscSqrtScalar(dp); /* beta = sqrt(dp); */
143: /* QR factorization */
144: coold = cold;
145: cold = c;
146: soold = sold;
147: sold = s;
148: rho0 = cold * alpha - coold * sold * betaold; /* gamma_bar */
149: rho1 = PetscSqrtScalar(rho0 * rho0 + beta * beta); /* gamma */
150: rho2 = sold * alpha + coold * cold * betaold; /* delta */
151: rho3 = soold * betaold; /* epsilon */
153: /* Givens rotation: [c -s; s c] (different from the Reference!) */
154: c = rho0 / rho1;
155: s = beta / rho1;
157: if (ksp->its == 1) ceta = beta1 / rho1;
158: else ceta = -(rho2 * ceta_old + rho3 * ceta_oold) / rho1;
160: s_prod = s_prod * PetscAbsScalar(s);
161: if (c == 0.0) np = s_prod * 1.e16;
162: else np = s_prod / PetscAbsScalar(c); /* residual norm for xc_k (CGNORM) */
164: if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = np;
165: else ksp->rnorm = 0.0;
166: KSPLogResidualHistory(ksp, ksp->rnorm);
167: KSPMonitor(ksp, i + 1, ksp->rnorm);
168: (*ksp->converged)(ksp, i + 1, ksp->rnorm, &ksp->reason, ksp->cnvP); /* test for convergence */
169: if (ksp->reason) break;
170: i++;
171: } while (i < ksp->max_it);
173: /* move to the CG point: xc_(k+1) */
174: if (c == 0.0) ceta_bar = ceta * 1.e15;
175: else ceta_bar = ceta / c;
177: VecAXPY(X, ceta_bar, Wbar); /* x <- x + ceta_bar*w_bar */
179: if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
180: return 0;
181: }
183: /*MC
184: KSPSYMMLQ - This code implements the SYMMLQ method.
186: Level: beginner
188: Notes:
189: The operator and the preconditioner must be symmetric for this method.
191: The preconditioner must be POSITIVE-DEFINITE.
193: Supports only left preconditioning.
195: Reference:
196: . * - Paige & Saunders, Solution of sparse indefinite systems of linear equations, SIAM J. Numer. Anal. 12, 1975.
198: .seealso: [](chapter_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`
199: M*/
200: PETSC_EXTERN PetscErrorCode KSPCreate_SYMMLQ(KSP ksp)
201: {
202: KSP_SYMMLQ *symmlq;
204: KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 3);
205: KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_LEFT, 1);
207: PetscNew(&symmlq);
208: symmlq->haptol = 1.e-18;
209: ksp->data = (void *)symmlq;
211: /*
212: Sets the functions that are associated with this data structure
213: (in C++ this is the same as defining virtual functions)
214: */
215: ksp->ops->setup = KSPSetUp_SYMMLQ;
216: ksp->ops->solve = KSPSolve_SYMMLQ;
217: ksp->ops->destroy = KSPDestroyDefault;
218: ksp->ops->setfromoptions = NULL;
219: ksp->ops->buildsolution = KSPBuildSolutionDefault;
220: ksp->ops->buildresidual = KSPBuildResidualDefault;
221: return 0;
222: }