Actual source code: qmrcgs.c
2: /*
3: This file implements QMRCGS (QMRCGStab).
4: Only right preconditioning is supported.
6: Contributed by: Xiangmin Jiao (xiangmin.jiao@stonybrook.edu)
8: References:
9: . * - Chan, Gallopoulos, Simoncini, Szeto, and Tong (SISC 1994), Ghai, Lu, and Jiao (NLAA 2019)
10: */
11: #include <../src/ksp/ksp/impls/bcgs/bcgsimpl.h>
13: static PetscErrorCode KSPSetUp_QMRCGS(KSP ksp)
14: {
15: KSPSetWorkVecs(ksp, 14);
16: return 0;
17: }
19: /* Only need a few hacks from KSPSolve_BCGS */
21: static PetscErrorCode KSPSolve_QMRCGS(KSP ksp)
22: {
23: PetscInt i;
24: PetscScalar eta, rho1, rho2, alpha, eta2, omega, beta, cf, cf1, uu;
25: Vec X, B, R, P, PH, V, D2, X2, S, SH, T, D, S2, RP, AX, Z;
26: PetscReal dp = 0.0, final, tau, tau2, theta, theta2, c, F, NV, vv;
27: KSP_BCGS *bcgs = (KSP_BCGS *)ksp->data;
28: PC pc;
29: Mat mat;
31: X = ksp->vec_sol;
32: B = ksp->vec_rhs;
33: R = ksp->work[0];
34: P = ksp->work[1];
35: PH = ksp->work[2];
36: V = ksp->work[3];
37: D2 = ksp->work[4];
38: X2 = ksp->work[5];
39: S = ksp->work[6];
40: SH = ksp->work[7];
41: T = ksp->work[8];
42: D = ksp->work[9];
43: S2 = ksp->work[10];
44: RP = ksp->work[11];
45: AX = ksp->work[12];
46: Z = ksp->work[13];
48: /* Only supports right preconditioning */
50: if (!ksp->guess_zero) {
51: if (!bcgs->guess) VecDuplicate(X, &bcgs->guess);
52: VecCopy(X, bcgs->guess);
53: } else {
54: VecSet(X, 0.0);
55: }
57: /* Compute initial residual */
58: KSPGetPC(ksp, &pc);
59: PCSetUp(pc);
60: PCGetOperators(pc, &mat, NULL);
61: if (!ksp->guess_zero) {
62: KSP_MatMult(ksp, mat, X, S2);
63: VecCopy(B, R);
64: VecAXPY(R, -1.0, S2);
65: } else {
66: VecCopy(B, R);
67: }
69: /* Test for nothing to do */
70: if (ksp->normtype != KSP_NORM_NONE) VecNorm(R, NORM_2, &dp);
71: PetscObjectSAWsTakeAccess((PetscObject)ksp);
72: ksp->its = 0;
73: ksp->rnorm = dp;
74: PetscObjectSAWsGrantAccess((PetscObject)ksp);
75: KSPLogResidualHistory(ksp, dp);
76: KSPMonitor(ksp, 0, dp);
77: (*ksp->converged)(ksp, 0, dp, &ksp->reason, ksp->cnvP);
78: if (ksp->reason) return 0;
80: /* Make the initial Rp == R */
81: VecCopy(R, RP);
83: eta = 1.0;
84: theta = 1.0;
85: if (dp == 0.0) {
86: VecNorm(R, NORM_2, &tau);
87: } else {
88: tau = dp;
89: }
91: VecDot(RP, RP, &rho1);
92: VecCopy(R, P);
94: i = 0;
95: do {
96: KSP_PCApply(ksp, P, PH); /* ph <- K p */
97: KSP_MatMult(ksp, mat, PH, V); /* v <- A ph */
99: VecDot(V, RP, &rho2); /* rho2 <- (v,rp) */
100: if (rho2 == 0.0) {
102: ksp->reason = KSP_DIVERGED_NANORINF;
103: break;
104: }
106: if (rho1 == 0) {
108: ksp->reason = KSP_DIVERGED_BREAKDOWN; /* Stagnation */
109: break;
110: }
112: alpha = rho1 / rho2;
113: VecWAXPY(S, -alpha, V, R); /* s <- r - alpha v */
115: /* First quasi-minimization step */
116: VecNorm(S, NORM_2, &F); /* f <- norm(s) */
117: theta2 = F / tau;
119: c = 1.0 / PetscSqrtReal(1.0 + theta2 * theta2);
121: tau2 = tau * theta2 * c;
122: eta2 = c * c * alpha;
123: cf = theta * theta * eta / alpha;
124: VecWAXPY(D2, cf, D, PH); /* d2 <- ph + cf d */
125: VecWAXPY(X2, eta2, D2, X); /* x2 <- x + eta2 d2 */
127: /* Apply the right preconditioner again */
128: KSP_PCApply(ksp, S, SH); /* sh <- K s */
129: KSP_MatMult(ksp, mat, SH, T); /* t <- A sh */
131: VecDotNorm2(S, T, &uu, &vv);
132: if (vv == 0.0) {
133: VecDot(S, S, &uu);
134: if (uu != 0.0) {
136: ksp->reason = KSP_DIVERGED_NANORINF;
137: break;
138: }
139: VecAXPY(X, alpha, SH); /* x <- x + alpha sh */
140: PetscObjectSAWsTakeAccess((PetscObject)ksp);
141: ksp->its++;
142: ksp->rnorm = 0.0;
143: ksp->reason = KSP_CONVERGED_RTOL;
144: PetscObjectSAWsGrantAccess((PetscObject)ksp);
145: KSPLogResidualHistory(ksp, dp);
146: KSPMonitor(ksp, i + 1, 0.0);
147: break;
148: }
149: VecNorm(V, NORM_2, &NV); /* nv <- norm(v) */
151: if (NV == 0) {
153: ksp->reason = KSP_DIVERGED_BREAKDOWN;
154: break;
155: }
157: if (uu == 0) {
159: ksp->reason = KSP_DIVERGED_BREAKDOWN; /* Stagnation */
160: break;
161: }
162: omega = uu / vv; /* omega <- uu/vv; */
164: /* Computing the residual */
165: VecWAXPY(R, -omega, T, S); /* r <- s - omega t */
167: /* Second quasi-minimization step */
168: if (ksp->normtype != KSP_NORM_NONE && ksp->chknorm < i + 2) VecNorm(R, NORM_2, &dp);
170: if (tau2 == 0) {
172: ksp->reason = KSP_DIVERGED_NANORINF;
173: break;
174: }
175: theta = dp / tau2;
176: c = 1.0 / PetscSqrtReal(1.0 + theta * theta);
177: if (dp == 0.0) {
178: VecNorm(R, NORM_2, &tau);
179: } else {
180: tau = dp;
181: }
182: tau = tau * c;
183: eta = c * c * omega;
185: cf1 = theta2 * theta2 * eta2 / omega;
186: VecWAXPY(D, cf1, D2, SH); /* d <- sh + cf1 d2 */
187: VecWAXPY(X, eta, D, X2); /* x <- x2 + eta d */
189: VecDot(R, RP, &rho2);
190: PetscObjectSAWsTakeAccess((PetscObject)ksp);
191: ksp->its++;
192: ksp->rnorm = dp;
193: PetscObjectSAWsGrantAccess((PetscObject)ksp);
194: KSPLogResidualHistory(ksp, dp);
195: KSPMonitor(ksp, i + 1, dp);
197: beta = (alpha * rho2) / (omega * rho1);
198: VecAXPBYPCZ(P, 1.0, -omega * beta, beta, R, V); /* p <- r - omega * beta* v + beta * p */
199: rho1 = rho2;
200: KSP_MatMult(ksp, mat, X, AX); /* Ax <- A x */
201: VecWAXPY(Z, -1.0, AX, B); /* r <- b - Ax */
202: VecNorm(Z, NORM_2, &final);
203: (*ksp->converged)(ksp, i + 1, dp, &ksp->reason, ksp->cnvP);
204: if (ksp->reason) break;
205: i++;
206: } while (i < ksp->max_it);
208: /* mark lack of convergence */
209: if (ksp->its >= ksp->max_it && !ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
210: return 0;
211: }
213: /*MC
214: KSPQMRCGS - Implements the QMRCGStab method.
216: Options Database Keys:
217: see KSPSolve()
219: Level: beginner
221: Notes:
222: Only right preconditioning is supported.
224: References:
225: . * - Chan, Gallopoulos, Simoncini, Szeto, and Tong (SISC 1994), Ghai, Lu, and Jiao (NLAA 2019)
227: .seealso: `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPBICG`, `KSPFBICGS`, `KSPFBCGSL`, `KSPSetPCSide()`
228: M*/
229: PETSC_EXTERN PetscErrorCode KSPCreate_QMRCGS(KSP ksp)
230: {
231: KSP_BCGS *bcgs;
232: static const char citations[] = "@article{chan1994qmrcgs,\n"
233: " title={A quasi-minimal residual variant of the {Bi-CGSTAB} algorithm for nonsymmetric systems},\n"
234: " author={Chan, Tony F and Gallopoulos, Efstratios and Simoncini, Valeria and Szeto, Tedd and Tong, Charles H},\n"
235: " journal={SIAM Journal on Scientific Computing},\n"
236: " volume={15},\n"
237: " number={2},\n"
238: " pages={338--347},\n"
239: " year={1994},\n"
240: " publisher={SIAM}\n"
241: "}\n"
242: "@article{ghai2019comparison,\n"
243: " title={A comparison of preconditioned {K}rylov subspace methods for large-scale nonsymmetric linear systems},\n"
244: " author={Ghai, Aditi and Lu, Cao and Jiao, Xiangmin},\n"
245: " journal={Numerical Linear Algebra with Applications},\n"
246: " volume={26},\n"
247: " number={1},\n"
248: " pages={e2215},\n"
249: " year={2019},\n"
250: " publisher={Wiley Online Library}\n"
251: "}\n";
252: PetscBool cite = PETSC_FALSE;
255: PetscCitationsRegister(citations, &cite);
256: PetscNew(&bcgs);
258: ksp->data = bcgs;
259: ksp->ops->setup = KSPSetUp_QMRCGS;
260: ksp->ops->solve = KSPSolve_QMRCGS;
261: ksp->ops->destroy = KSPDestroy_BCGS;
262: ksp->ops->reset = KSPReset_BCGS;
263: ksp->ops->buildresidual = KSPBuildResidualDefault;
264: ksp->ops->setfromoptions = KSPSetFromOptions_BCGS;
265: ksp->pc_side = PC_RIGHT; /* set default PC side */
267: KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_RIGHT, 2);
268: KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_RIGHT, 1);
269: return 0;
270: }