Actual source code: bcgs.c
2: #include <../src/ksp/ksp/impls/bcgs/bcgsimpl.h>
4: PetscErrorCode KSPSetFromOptions_BCGS(KSP ksp, PetscOptionItems *PetscOptionsObject)
5: {
6: PetscOptionsHeadBegin(PetscOptionsObject, "KSP BCGS Options");
7: PetscOptionsHeadEnd();
8: return 0;
9: }
11: PetscErrorCode KSPSetUp_BCGS(KSP ksp)
12: {
13: KSPSetWorkVecs(ksp, 6);
14: return 0;
15: }
17: PetscErrorCode KSPSolve_BCGS(KSP ksp)
18: {
19: PetscInt i;
20: PetscScalar rho, rhoold, alpha, beta, omega, omegaold, d1;
21: Vec X, B, V, P, R, RP, T, S;
22: PetscReal dp = 0.0, d2;
23: KSP_BCGS *bcgs = (KSP_BCGS *)ksp->data;
25: X = ksp->vec_sol;
26: B = ksp->vec_rhs;
27: R = ksp->work[0];
28: RP = ksp->work[1];
29: V = ksp->work[2];
30: T = ksp->work[3];
31: S = ksp->work[4];
32: P = ksp->work[5];
34: /* Compute initial preconditioned residual */
35: KSPInitialResidual(ksp, X, V, T, R, B);
37: /* with right preconditioning need to save initial guess to add to final solution */
38: if (ksp->pc_side == PC_RIGHT && !ksp->guess_zero) {
39: if (!bcgs->guess) VecDuplicate(X, &bcgs->guess);
40: VecCopy(X, bcgs->guess);
41: VecSet(X, 0.0);
42: }
44: /* Test for nothing to do */
45: if (ksp->normtype != KSP_NORM_NONE) {
46: VecNorm(R, NORM_2, &dp);
47: KSPCheckNorm(ksp, dp);
48: }
49: PetscObjectSAWsTakeAccess((PetscObject)ksp);
50: ksp->its = 0;
51: ksp->rnorm = dp;
52: PetscObjectSAWsGrantAccess((PetscObject)ksp);
53: KSPLogResidualHistory(ksp, dp);
54: KSPMonitor(ksp, 0, dp);
55: (*ksp->converged)(ksp, 0, dp, &ksp->reason, ksp->cnvP);
56: if (ksp->reason) {
57: if (bcgs->guess) VecAXPY(X, 1.0, bcgs->guess);
58: return 0;
59: }
61: /* Make the initial Rp == R */
62: VecCopy(R, RP);
64: rhoold = 1.0;
65: alpha = 1.0;
66: omegaold = 1.0;
67: VecSet(P, 0.0);
68: VecSet(V, 0.0);
70: i = 0;
71: do {
72: VecDot(R, RP, &rho); /* rho <- (r,rp) */
73: beta = (rho / rhoold) * (alpha / omegaold);
74: VecAXPBYPCZ(P, 1.0, -omegaold * beta, beta, R, V); /* p <- r - omega * beta* v + beta * p */
75: KSP_PCApplyBAorAB(ksp, P, V, T); /* v <- K p */
76: VecDot(V, RP, &d1);
77: KSPCheckDot(ksp, d1);
78: if (d1 == 0.0) {
80: ksp->reason = KSP_DIVERGED_BREAKDOWN;
81: PetscInfo(ksp, "Breakdown due to zero inner product\n");
82: break;
83: }
84: alpha = rho / d1; /* a <- rho / (v,rp) */
85: VecWAXPY(S, -alpha, V, R); /* s <- r - a v */
86: KSP_PCApplyBAorAB(ksp, S, T, R); /* t <- K s */
87: VecDotNorm2(S, T, &d1, &d2);
88: if (d2 == 0.0) {
89: /* t is 0. if s is 0, then alpha v == r, and hence alpha p
90: may be our solution. Give it a try? */
91: VecDot(S, S, &d1);
92: if (d1 != 0.0) {
94: ksp->reason = KSP_DIVERGED_BREAKDOWN;
95: PetscInfo(ksp, "Failed due to singular preconditioned operator\n");
96: break;
97: }
98: VecAXPY(X, alpha, P); /* x <- x + a p */
99: PetscObjectSAWsTakeAccess((PetscObject)ksp);
100: ksp->its++;
101: ksp->rnorm = 0.0;
102: ksp->reason = KSP_CONVERGED_RTOL;
103: PetscObjectSAWsGrantAccess((PetscObject)ksp);
104: KSPLogResidualHistory(ksp, dp);
105: KSPMonitor(ksp, i + 1, 0.0);
106: break;
107: }
108: omega = d1 / d2; /* w <- (t's) / (t't) */
109: VecAXPBYPCZ(X, alpha, omega, 1.0, P, S); /* x <- alpha * p + omega * s + x */
110: VecWAXPY(R, -omega, T, S); /* r <- s - w t */
111: if (ksp->normtype != KSP_NORM_NONE && ksp->chknorm < i + 2) {
112: VecNorm(R, NORM_2, &dp);
113: KSPCheckNorm(ksp, dp);
114: }
116: rhoold = rho;
117: omegaold = omega;
119: PetscObjectSAWsTakeAccess((PetscObject)ksp);
120: ksp->its++;
121: ksp->rnorm = dp;
122: PetscObjectSAWsGrantAccess((PetscObject)ksp);
123: KSPLogResidualHistory(ksp, dp);
124: KSPMonitor(ksp, i + 1, dp);
125: (*ksp->converged)(ksp, i + 1, dp, &ksp->reason, ksp->cnvP);
126: if (ksp->reason) break;
127: if (rho == 0.0) {
129: ksp->reason = KSP_DIVERGED_BREAKDOWN;
130: PetscInfo(ksp, "Breakdown due to zero rho inner product\n");
131: break;
132: }
133: i++;
134: } while (i < ksp->max_it);
136: if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
138: KSPUnwindPreconditioner(ksp, X, T);
139: if (bcgs->guess) VecAXPY(X, 1.0, bcgs->guess);
140: return 0;
141: }
143: PetscErrorCode KSPBuildSolution_BCGS(KSP ksp, Vec v, Vec *V)
144: {
145: KSP_BCGS *bcgs = (KSP_BCGS *)ksp->data;
147: if (ksp->pc_side == PC_RIGHT) {
148: if (v) {
149: KSP_PCApply(ksp, ksp->vec_sol, v);
150: if (bcgs->guess) VecAXPY(v, 1.0, bcgs->guess);
151: *V = v;
152: } else SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Not working with right preconditioner");
153: } else {
154: if (v) {
155: VecCopy(ksp->vec_sol, v);
156: *V = v;
157: } else *V = ksp->vec_sol;
158: }
159: return 0;
160: }
162: PetscErrorCode KSPReset_BCGS(KSP ksp)
163: {
164: KSP_BCGS *cg = (KSP_BCGS *)ksp->data;
166: VecDestroy(&cg->guess);
167: return 0;
168: }
170: PetscErrorCode KSPDestroy_BCGS(KSP ksp)
171: {
172: KSPReset_BCGS(ksp);
173: KSPDestroyDefault(ksp);
174: return 0;
175: }
177: /*MC
178: KSPBCGS - Implements the BiCGStab (Stabilized version of Biconjugate Gradient) method.
180: Level: beginner
182: Notes:
183: Supports left and right preconditioning but not symmetric
185: See `KSPBCGSL` for additional stabilization
187: See `KSPFBCGS`, `KSPFBCGSR`, and `KSPPIPEBCGS` for flexible and pipelined versions of the algorithm
189: Reference:
190: . * - van der Vorst, SIAM J. Sci. Stat. Comput., 1992.
192: .seealso: [](chapter_ksp), `KSPFBCGS`, `KSPFBCGSR`, `KSPPIPEBCGS`, `KSPBCGSL`, `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPBICG`, `KSPBCGSL`, `KSPFBICG`, `KSPQMRCGS`, `KSPSetPCSide()`
193: M*/
194: PETSC_EXTERN PetscErrorCode KSPCreate_BCGS(KSP ksp)
195: {
196: KSP_BCGS *bcgs;
198: PetscNew(&bcgs);
200: ksp->data = bcgs;
201: ksp->ops->setup = KSPSetUp_BCGS;
202: ksp->ops->solve = KSPSolve_BCGS;
203: ksp->ops->destroy = KSPDestroy_BCGS;
204: ksp->ops->reset = KSPReset_BCGS;
205: ksp->ops->buildsolution = KSPBuildSolution_BCGS;
206: ksp->ops->buildresidual = KSPBuildResidualDefault;
207: ksp->ops->setfromoptions = KSPSetFromOptions_BCGS;
209: KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 3);
210: KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_RIGHT, 2);
211: KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_LEFT, 1);
212: KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_RIGHT, 1);
213: return 0;
214: }