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