Actual source code: fbcgs.c
2: /*
3: This file implements flexible BiCGStab (FBiCGStab).
4: Only allow right preconditioning.
5: */
6: #include <../src/ksp/ksp/impls/bcgs/bcgsimpl.h>
8: static PetscErrorCode KSPSetUp_FBCGS(KSP ksp)
9: {
10: KSPSetWorkVecs(ksp, 8);
11: return 0;
12: }
14: /* Only need a few hacks from KSPSolve_BCGS */
16: static PetscErrorCode KSPSolve_FBCGS(KSP ksp)
17: {
18: PetscInt i;
19: PetscScalar rho, rhoold, alpha, beta, omega, omegaold, d1;
20: Vec X, B, V, P, R, RP, T, S, P2, S2;
21: PetscReal dp = 0.0, d2;
22: KSP_BCGS *bcgs = (KSP_BCGS *)ksp->data;
23: PC pc;
24: Mat mat;
26: X = ksp->vec_sol;
27: B = ksp->vec_rhs;
28: R = ksp->work[0];
29: RP = ksp->work[1];
30: V = ksp->work[2];
31: T = ksp->work[3];
32: S = ksp->work[4];
33: P = ksp->work[5];
34: S2 = ksp->work[6];
35: P2 = ksp->work[7];
37: /* Only supports right preconditioning */
39: if (!ksp->guess_zero) {
40: if (!bcgs->guess) VecDuplicate(X, &bcgs->guess);
41: VecCopy(X, bcgs->guess);
42: } else {
43: VecSet(X, 0.0);
44: }
46: /* Compute initial residual */
47: KSPGetPC(ksp, &pc);
48: PCSetUp(pc);
49: PCGetOperators(pc, &mat, NULL);
50: if (!ksp->guess_zero) {
51: KSP_MatMult(ksp, mat, X, S2);
52: VecCopy(B, R);
53: VecAXPY(R, -1.0, S2);
54: } else {
55: VecCopy(B, R);
56: }
58: /* Test for nothing to do */
59: if (ksp->normtype != KSP_NORM_NONE) VecNorm(R, NORM_2, &dp);
60: PetscObjectSAWsTakeAccess((PetscObject)ksp);
61: ksp->its = 0;
62: ksp->rnorm = dp;
63: PetscObjectSAWsGrantAccess((PetscObject)ksp);
64: KSPLogResidualHistory(ksp, dp);
65: KSPMonitor(ksp, 0, dp);
66: (*ksp->converged)(ksp, 0, dp, &ksp->reason, ksp->cnvP);
67: if (ksp->reason) return 0;
69: /* Make the initial Rp == R */
70: VecCopy(R, RP);
72: rhoold = 1.0;
73: alpha = 1.0;
74: omegaold = 1.0;
75: VecSet(P, 0.0);
76: VecSet(V, 0.0);
78: i = 0;
79: do {
80: VecDot(R, RP, &rho); /* rho <- (r,rp) */
81: beta = (rho / rhoold) * (alpha / omegaold);
82: VecAXPBYPCZ(P, 1.0, -omegaold * beta, beta, R, V); /* p <- r - omega * beta* v + beta * p */
84: KSP_PCApply(ksp, P, P2); /* p2 <- K p */
85: KSP_MatMult(ksp, mat, P2, V); /* v <- A p2 */
87: VecDot(V, RP, &d1);
88: if (d1 == 0.0) {
90: ksp->reason = KSP_DIVERGED_BREAKDOWN;
91: PetscInfo(ksp, "Breakdown due to zero inner product\n");
92: break;
93: }
94: alpha = rho / d1; /* alpha <- rho / (v,rp) */
95: VecWAXPY(S, -alpha, V, R); /* s <- r - alpha v */
97: KSP_PCApply(ksp, S, S2); /* s2 <- K s */
98: KSP_MatMult(ksp, mat, S2, T); /* t <- A s2 */
100: VecDotNorm2(S, T, &d1, &d2);
101: if (d2 == 0.0) {
102: /* t is 0. if s is 0, then alpha v == r, and hence alpha p may be our solution. Give it a try? */
103: VecDot(S, S, &d1);
104: if (d1 != 0.0) {
106: ksp->reason = KSP_DIVERGED_BREAKDOWN;
107: PetscInfo(ksp, "Failed due to singular preconditioned operator\n");
108: break;
109: }
110: VecAXPY(X, alpha, P2); /* x <- x + alpha p2 */
111: PetscObjectSAWsTakeAccess((PetscObject)ksp);
112: ksp->its++;
113: ksp->rnorm = 0.0;
114: ksp->reason = KSP_CONVERGED_RTOL;
115: PetscObjectSAWsGrantAccess((PetscObject)ksp);
116: KSPLogResidualHistory(ksp, dp);
117: KSPMonitor(ksp, i + 1, 0.0);
118: break;
119: }
120: omega = d1 / d2; /* omega <- (t's) / (t't) */
121: VecAXPBYPCZ(X, alpha, omega, 1.0, P2, S2); /* x <- alpha * p2 + omega * s2 + x */
123: VecWAXPY(R, -omega, T, S); /* r <- s - omega t */
124: if (ksp->normtype != KSP_NORM_NONE && ksp->chknorm < i + 2) VecNorm(R, NORM_2, &dp);
126: rhoold = rho;
127: omegaold = omega;
129: PetscObjectSAWsTakeAccess((PetscObject)ksp);
130: ksp->its++;
131: ksp->rnorm = dp;
132: PetscObjectSAWsGrantAccess((PetscObject)ksp);
133: KSPLogResidualHistory(ksp, dp);
134: KSPMonitor(ksp, i + 1, dp);
135: (*ksp->converged)(ksp, i + 1, dp, &ksp->reason, ksp->cnvP);
136: if (ksp->reason) break;
137: if (rho == 0.0) {
139: ksp->reason = KSP_DIVERGED_BREAKDOWN;
140: PetscInfo(ksp, "Breakdown due to zero rho inner product\n");
141: break;
142: }
143: i++;
144: } while (i < ksp->max_it);
146: if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
147: return 0;
148: }
150: /*MC
151: KSPFBCGS - Implements flexible BiCGStab method.
153: Options Database Keys:
154: see KSPSolve()
156: Level: beginner
158: Notes:
159: Only allow right preconditioning
161: .seealso: `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPBICG`, `KSPFBCGSL`, `KSPSetPCSide()`
162: M*/
163: PETSC_EXTERN PetscErrorCode KSPCreate_FBCGS(KSP ksp)
164: {
165: KSP_BCGS *bcgs;
167: PetscNew(&bcgs);
169: ksp->data = bcgs;
170: ksp->ops->setup = KSPSetUp_FBCGS;
171: ksp->ops->solve = KSPSolve_FBCGS;
172: ksp->ops->destroy = KSPDestroy_BCGS;
173: ksp->ops->reset = KSPReset_BCGS;
174: ksp->ops->buildresidual = KSPBuildResidualDefault;
175: ksp->ops->setfromoptions = KSPSetFromOptions_BCGS;
176: ksp->pc_side = PC_RIGHT; /* set default PC side */
178: KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 3);
179: KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_RIGHT, 2);
180: KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_LEFT, 1);
181: KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_RIGHT, 1);
182: return 0;
183: }