Actual source code: bcgsl.c
1: /*
2: * Implementation of BiCGstab(L) the paper by D.R. Fokkema,
3: * "Enhanced implementation of BiCGStab(L) for solving linear systems
4: * of equations". This uses tricky delayed updating ideas to prevent
5: * round-off buildup.
6: *
7: * This has not been completely cleaned up into PETSc style.
8: *
9: * All the BLAS and LAPACK calls below should be removed and replaced with
10: * loops and the macros for block solvers converted from LINPACK; there is no way
11: * calls to BLAS/LAPACK make sense for size 2, 3, 4, etc.
12: */
13: #include <petsc/private/kspimpl.h>
14: #include <../src/ksp/ksp/impls/bcgsl/bcgslimpl.h>
15: #include <petscblaslapack.h>
17: static PetscErrorCode KSPSolve_BCGSL(KSP ksp)
18: {
19: KSP_BCGSL *bcgsl = (KSP_BCGSL *)ksp->data;
20: PetscScalar alpha, beta, omega, sigma;
21: PetscScalar rho0, rho1;
22: PetscReal kappa0, kappaA, kappa1;
23: PetscReal ghat;
24: PetscReal zeta, zeta0, rnmax_computed, rnmax_true, nrm0;
25: PetscBool bUpdateX;
26: PetscInt maxit;
27: PetscInt h, i, j, k, vi, ell;
28: PetscBLASInt ldMZ, bierr;
29: PetscScalar utb;
30: PetscReal max_s, pinv_tol;
32: /* set up temporary vectors */
33: vi = 0;
34: ell = bcgsl->ell;
35: bcgsl->vB = ksp->work[vi];
36: vi++;
37: bcgsl->vRt = ksp->work[vi];
38: vi++;
39: bcgsl->vTm = ksp->work[vi];
40: vi++;
41: bcgsl->vvR = ksp->work + vi;
42: vi += ell + 1;
43: bcgsl->vvU = ksp->work + vi;
44: vi += ell + 1;
45: bcgsl->vXr = ksp->work[vi];
46: vi++;
47: PetscBLASIntCast(ell + 1, &ldMZ);
49: /* Prime the iterative solver */
50: KSPInitialResidual(ksp, VX, VTM, VB, VVR[0], ksp->vec_rhs);
51: VecNorm(VVR[0], NORM_2, &zeta0);
52: KSPCheckNorm(ksp, zeta0);
53: rnmax_computed = zeta0;
54: rnmax_true = zeta0;
56: PetscObjectSAWsTakeAccess((PetscObject)ksp);
57: ksp->its = 0;
58: if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = zeta0;
59: else ksp->rnorm = 0.0;
60: PetscObjectSAWsGrantAccess((PetscObject)ksp);
61: (*ksp->converged)(ksp, 0, ksp->rnorm, &ksp->reason, ksp->cnvP);
62: if (ksp->reason) return 0;
64: VecSet(VVU[0], 0.0);
65: alpha = 0.;
66: rho0 = omega = 1;
68: if (bcgsl->delta > 0.0) {
69: VecCopy(VX, VXR);
70: VecSet(VX, 0.0);
71: VecCopy(VVR[0], VB);
72: } else {
73: VecCopy(ksp->vec_rhs, VB);
74: }
76: /* Life goes on */
77: VecCopy(VVR[0], VRT);
78: zeta = zeta0;
80: KSPGetTolerances(ksp, NULL, NULL, NULL, &maxit);
82: for (k = 0; k < maxit; k += bcgsl->ell) {
83: ksp->its = k;
84: if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = zeta;
85: else ksp->rnorm = 0.0;
87: KSPLogResidualHistory(ksp, ksp->rnorm);
88: KSPMonitor(ksp, ksp->its, ksp->rnorm);
90: (*ksp->converged)(ksp, k, ksp->rnorm, &ksp->reason, ksp->cnvP);
91: if (ksp->reason < 0) return 0;
92: if (ksp->reason) {
93: if (bcgsl->delta > 0.0) VecAXPY(VX, 1.0, VXR);
94: return 0;
95: }
97: /* BiCG part */
98: rho0 = -omega * rho0;
99: nrm0 = zeta;
100: for (j = 0; j < bcgsl->ell; j++) {
101: /* rho1 <- r_j' * r_tilde */
102: VecDot(VVR[j], VRT, &rho1);
103: KSPCheckDot(ksp, rho1);
104: if (rho1 == 0.0) {
105: ksp->reason = KSP_DIVERGED_BREAKDOWN_BICG;
106: return 0;
107: }
108: beta = alpha * (rho1 / rho0);
109: rho0 = rho1;
110: for (i = 0; i <= j; i++) {
111: /* u_i <- r_i - beta*u_i */
112: VecAYPX(VVU[i], -beta, VVR[i]);
113: }
114: /* u_{j+1} <- inv(K)*A*u_j */
115: KSP_PCApplyBAorAB(ksp, VVU[j], VVU[j + 1], VTM);
117: VecDot(VVU[j + 1], VRT, &sigma);
118: KSPCheckDot(ksp, sigma);
119: if (sigma == 0.0) {
120: ksp->reason = KSP_DIVERGED_BREAKDOWN_BICG;
121: return 0;
122: }
123: alpha = rho1 / sigma;
125: /* x <- x + alpha*u_0 */
126: VecAXPY(VX, alpha, VVU[0]);
128: for (i = 0; i <= j; i++) {
129: /* r_i <- r_i - alpha*u_{i+1} */
130: VecAXPY(VVR[i], -alpha, VVU[i + 1]);
131: }
133: /* r_{j+1} <- inv(K)*A*r_j */
134: KSP_PCApplyBAorAB(ksp, VVR[j], VVR[j + 1], VTM);
136: VecNorm(VVR[0], NORM_2, &nrm0);
137: KSPCheckNorm(ksp, nrm0);
138: if (bcgsl->delta > 0.0) {
139: if (rnmax_computed < nrm0) rnmax_computed = nrm0;
140: if (rnmax_true < nrm0) rnmax_true = nrm0;
141: }
143: /* NEW: check for early exit */
144: (*ksp->converged)(ksp, k + j, nrm0, &ksp->reason, ksp->cnvP);
145: if (ksp->reason) {
146: PetscObjectSAWsTakeAccess((PetscObject)ksp);
147: ksp->its = k + j;
148: ksp->rnorm = nrm0;
150: PetscObjectSAWsGrantAccess((PetscObject)ksp);
151: if (ksp->reason < 0) return 0;
152: }
153: }
155: /* Polynomial part */
156: for (i = 0; i <= bcgsl->ell; ++i) VecMDot(VVR[i], i + 1, VVR, &MZa[i * ldMZ]);
157: /* Symmetrize MZa */
158: for (i = 0; i <= bcgsl->ell; ++i) {
159: for (j = i + 1; j <= bcgsl->ell; ++j) MZa[i * ldMZ + j] = MZa[j * ldMZ + i] = PetscConj(MZa[j * ldMZ + i]);
160: }
161: /* Copy MZa to MZb */
162: PetscArraycpy(MZb, MZa, ldMZ * ldMZ);
164: if (!bcgsl->bConvex || bcgsl->ell == 1) {
165: PetscBLASInt ione = 1, bell;
166: PetscBLASIntCast(bcgsl->ell, &bell);
168: AY0c[0] = -1;
169: if (bcgsl->pinv) {
170: #if defined(PETSC_USE_COMPLEX)
171: PetscCallBLAS("LAPACKgesvd", LAPACKgesvd_("A", "A", &bell, &bell, &MZa[1 + ldMZ], &ldMZ, bcgsl->s, bcgsl->u, &bell, bcgsl->v, &bell, bcgsl->work, &bcgsl->lwork, bcgsl->realwork, &bierr));
172: #else
173: PetscCallBLAS("LAPACKgesvd", LAPACKgesvd_("A", "A", &bell, &bell, &MZa[1 + ldMZ], &ldMZ, bcgsl->s, bcgsl->u, &bell, bcgsl->v, &bell, bcgsl->work, &bcgsl->lwork, &bierr));
174: #endif
175: if (bierr != 0) {
176: ksp->reason = KSP_DIVERGED_BREAKDOWN;
177: return 0;
178: }
179: /* Apply pseudo-inverse */
180: max_s = bcgsl->s[0];
181: for (i = 1; i < bell; i++) {
182: if (bcgsl->s[i] > max_s) max_s = bcgsl->s[i];
183: }
184: /* tolerance is hardwired to bell*max(s)*PETSC_MACHINE_EPSILON */
185: pinv_tol = bell * max_s * PETSC_MACHINE_EPSILON;
186: PetscArrayzero(&AY0c[1], bell);
187: for (i = 0; i < bell; i++) {
188: if (bcgsl->s[i] >= pinv_tol) {
189: utb = 0.;
190: for (j = 0; j < bell; j++) utb += MZb[1 + j] * bcgsl->u[i * bell + j];
192: for (j = 0; j < bell; j++) AY0c[1 + j] += utb / bcgsl->s[i] * bcgsl->v[j * bell + i];
193: }
194: }
195: } else {
196: PetscCallBLAS("LAPACKpotrf", LAPACKpotrf_("Lower", &bell, &MZa[1 + ldMZ], &ldMZ, &bierr));
197: if (bierr != 0) {
198: ksp->reason = KSP_DIVERGED_BREAKDOWN;
199: return 0;
200: }
201: PetscArraycpy(&AY0c[1], &MZb[1], bcgsl->ell);
202: PetscCallBLAS("LAPACKpotrs", LAPACKpotrs_("Lower", &bell, &ione, &MZa[1 + ldMZ], &ldMZ, &AY0c[1], &ldMZ, &bierr));
203: }
204: } else {
205: PetscBLASInt ione = 1;
206: PetscScalar aone = 1.0, azero = 0.0;
207: PetscBLASInt neqs;
208: PetscBLASIntCast(bcgsl->ell - 1, &neqs);
210: PetscCallBLAS("LAPACKpotrf", LAPACKpotrf_("Lower", &neqs, &MZa[1 + ldMZ], &ldMZ, &bierr));
211: if (bierr != 0) {
212: ksp->reason = KSP_DIVERGED_BREAKDOWN;
213: return 0;
214: }
215: PetscArraycpy(&AY0c[1], &MZb[1], bcgsl->ell - 1);
216: PetscCallBLAS("LAPACKpotrs", LAPACKpotrs_("Lower", &neqs, &ione, &MZa[1 + ldMZ], &ldMZ, &AY0c[1], &ldMZ, &bierr));
217: AY0c[0] = -1;
218: AY0c[bcgsl->ell] = 0.;
220: PetscArraycpy(&AYlc[1], &MZb[1 + ldMZ * (bcgsl->ell)], bcgsl->ell - 1);
221: PetscCallBLAS("LAPACKpotrs", LAPACKpotrs_("Lower", &neqs, &ione, &MZa[1 + ldMZ], &ldMZ, &AYlc[1], &ldMZ, &bierr));
223: AYlc[0] = 0.;
224: AYlc[bcgsl->ell] = -1;
226: PetscCallBLAS("BLASgemv", BLASgemv_("NoTr", &ldMZ, &ldMZ, &aone, MZb, &ldMZ, AY0c, &ione, &azero, AYtc, &ione));
228: kappa0 = PetscRealPart(BLASdot_(&ldMZ, AY0c, &ione, AYtc, &ione));
230: /* round-off can cause negative kappa's */
231: if (kappa0 < 0) kappa0 = -kappa0;
232: kappa0 = PetscSqrtReal(kappa0);
234: kappaA = PetscRealPart(BLASdot_(&ldMZ, AYlc, &ione, AYtc, &ione));
236: PetscCallBLAS("BLASgemv", BLASgemv_("noTr", &ldMZ, &ldMZ, &aone, MZb, &ldMZ, AYlc, &ione, &azero, AYtc, &ione));
238: kappa1 = PetscRealPart(BLASdot_(&ldMZ, AYlc, &ione, AYtc, &ione));
240: if (kappa1 < 0) kappa1 = -kappa1;
241: kappa1 = PetscSqrtReal(kappa1);
243: if (kappa0 != 0.0 && kappa1 != 0.0) {
244: if (kappaA < 0.7 * kappa0 * kappa1) {
245: ghat = (kappaA < 0.0) ? -0.7 * kappa0 / kappa1 : 0.7 * kappa0 / kappa1;
246: } else {
247: ghat = kappaA / (kappa1 * kappa1);
248: }
249: for (i = 0; i <= bcgsl->ell; i++) AY0c[i] = AY0c[i] - ghat * AYlc[i];
250: }
251: }
253: omega = AY0c[bcgsl->ell];
254: for (h = bcgsl->ell; h > 0 && omega == 0.0; h--) omega = AY0c[h];
255: if (omega == 0.0) {
256: ksp->reason = KSP_DIVERGED_BREAKDOWN;
257: return 0;
258: }
260: VecMAXPY(VX, bcgsl->ell, AY0c + 1, VVR);
261: for (i = 1; i <= bcgsl->ell; i++) AY0c[i] *= -1.0;
262: VecMAXPY(VVU[0], bcgsl->ell, AY0c + 1, VVU + 1);
263: VecMAXPY(VVR[0], bcgsl->ell, AY0c + 1, VVR + 1);
264: for (i = 1; i <= bcgsl->ell; i++) AY0c[i] *= -1.0;
265: VecNorm(VVR[0], NORM_2, &zeta);
266: KSPCheckNorm(ksp, zeta);
268: /* Accurate Update */
269: if (bcgsl->delta > 0.0) {
270: if (rnmax_computed < zeta) rnmax_computed = zeta;
271: if (rnmax_true < zeta) rnmax_true = zeta;
273: bUpdateX = (PetscBool)(zeta < bcgsl->delta * zeta0 && zeta0 <= rnmax_computed);
274: if ((zeta < bcgsl->delta * rnmax_true && zeta0 <= rnmax_true) || bUpdateX) {
275: /* r0 <- b-inv(K)*A*X */
276: KSP_PCApplyBAorAB(ksp, VX, VVR[0], VTM);
277: VecAYPX(VVR[0], -1.0, VB);
278: rnmax_true = zeta;
280: if (bUpdateX) {
281: VecAXPY(VXR, 1.0, VX);
282: VecSet(VX, 0.0);
283: VecCopy(VVR[0], VB);
284: rnmax_computed = zeta;
285: }
286: }
287: }
288: }
289: if (bcgsl->delta > 0.0) VecAXPY(VX, 1.0, VXR);
291: ksp->its = k;
292: if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = zeta;
293: else ksp->rnorm = 0.0;
294: KSPMonitor(ksp, ksp->its, ksp->rnorm);
295: KSPLogResidualHistory(ksp, ksp->rnorm);
296: (*ksp->converged)(ksp, k, ksp->rnorm, &ksp->reason, ksp->cnvP);
297: if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
298: return 0;
299: }
301: /*@
302: KSPBCGSLSetXRes - Sets the parameter governing when
303: exact residuals will be used instead of computed residuals.
305: Logically Collective on ksp
307: Input Parameters:
308: + ksp - iterative context obtained from KSPCreate
309: - delta - computed residuals are used alone when delta is not positive
311: Options Database Keys:
313: . -ksp_bcgsl_xres delta - Threshold used to decide when to refresh computed residuals
315: Level: intermediate
317: .seealso: `KSPBCGSLSetEll()`, `KSPBCGSLSetPol()`, `KSP`
318: @*/
319: PetscErrorCode KSPBCGSLSetXRes(KSP ksp, PetscReal delta)
320: {
321: KSP_BCGSL *bcgsl = (KSP_BCGSL *)ksp->data;
324: if (ksp->setupstage) {
325: if ((delta <= 0 && bcgsl->delta > 0) || (delta > 0 && bcgsl->delta <= 0)) {
326: VecDestroyVecs(ksp->nwork, &ksp->work);
327: PetscFree5(AY0c, AYlc, AYtc, MZa, MZb);
328: PetscFree4(bcgsl->work, bcgsl->s, bcgsl->u, bcgsl->v);
329: ksp->setupstage = KSP_SETUP_NEW;
330: }
331: }
332: bcgsl->delta = delta;
333: return 0;
334: }
336: /*@
337: KSPBCGSLSetUsePseudoinverse - Use pseudoinverse (via SVD) to solve polynomial part of update
339: Logically Collective on ksp
341: Input Parameters:
342: + ksp - iterative context obtained from KSPCreate
343: - use_pinv - set to PETSC_TRUE when using pseudoinverse
345: Options Database Keys:
347: . -ksp_bcgsl_pinv - use pseudoinverse
349: Level: intermediate
351: .seealso: `KSPBCGSLSetEll()`, `KSP`
352: @*/
353: PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP ksp, PetscBool use_pinv)
354: {
355: KSP_BCGSL *bcgsl = (KSP_BCGSL *)ksp->data;
357: bcgsl->pinv = use_pinv;
358: return 0;
359: }
361: /*@
362: KSPBCGSLSetPol - Sets the type of polynomial part will
363: be used in the BiCGSTab(L) solver.
365: Logically Collective on ksp
367: Input Parameters:
368: + ksp - iterative context obtained from KSPCreate
369: - uMROR - set to PETSC_TRUE when the polynomial is a convex combination of an MR and an OR step.
371: Options Database Keys:
373: + -ksp_bcgsl_cxpoly - use enhanced polynomial
374: - -ksp_bcgsl_mrpoly - use standard polynomial
376: Level: intermediate
378: .seealso: `KSP`, `KSPBCGSL`, `KSPCreate()`, `KSPSetType()`
379: @*/
380: PetscErrorCode KSPBCGSLSetPol(KSP ksp, PetscBool uMROR)
381: {
382: KSP_BCGSL *bcgsl = (KSP_BCGSL *)ksp->data;
386: if (!ksp->setupstage) bcgsl->bConvex = uMROR;
387: else if (bcgsl->bConvex != uMROR) {
388: /* free the data structures,
389: then create them again
390: */
391: VecDestroyVecs(ksp->nwork, &ksp->work);
392: PetscFree5(AY0c, AYlc, AYtc, MZa, MZb);
393: PetscFree4(bcgsl->work, bcgsl->s, bcgsl->u, bcgsl->v);
395: bcgsl->bConvex = uMROR;
396: ksp->setupstage = KSP_SETUP_NEW;
397: }
398: return 0;
399: }
401: /*@
402: KSPBCGSLSetEll - Sets the number of search directions in BiCGStab(L).
404: Logically Collective on ksp
406: Input Parameters:
407: + ksp - iterative context obtained from KSPCreate
408: - ell - number of search directions
410: Options Database Keys:
412: . -ksp_bcgsl_ell ell - Number of Krylov search directions
414: Level: intermediate
416: Notes:
417: For large ell it is common for the polynomial update problem to become singular (due to happy breakdown for smallish
418: test problems, but also for larger problems). Consequently, by default, the system is solved by pseudoinverse, which
419: allows the iteration to complete successfully. See KSPBCGSLSetUsePseudoinverse() to switch to a conventional solve.
421: .seealso: `KSPBCGSLSetUsePseudoinverse()`, `KSP`, `KSPBCGSL`
422: @*/
423: PetscErrorCode KSPBCGSLSetEll(KSP ksp, PetscInt ell)
424: {
425: KSP_BCGSL *bcgsl = (KSP_BCGSL *)ksp->data;
430: if (!ksp->setupstage) bcgsl->ell = ell;
431: else if (bcgsl->ell != ell) {
432: /* free the data structures, then create them again */
433: VecDestroyVecs(ksp->nwork, &ksp->work);
434: PetscFree5(AY0c, AYlc, AYtc, MZa, MZb);
435: PetscFree4(bcgsl->work, bcgsl->s, bcgsl->u, bcgsl->v);
437: bcgsl->ell = ell;
438: ksp->setupstage = KSP_SETUP_NEW;
439: }
440: return 0;
441: }
443: PetscErrorCode KSPView_BCGSL(KSP ksp, PetscViewer viewer)
444: {
445: KSP_BCGSL *bcgsl = (KSP_BCGSL *)ksp->data;
446: PetscBool isascii;
448: PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);
450: if (isascii) {
451: PetscViewerASCIIPrintf(viewer, " Ell = %" PetscInt_FMT "\n", bcgsl->ell);
452: PetscViewerASCIIPrintf(viewer, " Delta = %g\n", (double)bcgsl->delta);
453: }
454: return 0;
455: }
457: PetscErrorCode KSPSetFromOptions_BCGSL(KSP ksp, PetscOptionItems *PetscOptionsObject)
458: {
459: KSP_BCGSL *bcgsl = (KSP_BCGSL *)ksp->data;
460: PetscInt this_ell;
461: PetscReal delta;
462: PetscBool flga = PETSC_FALSE, flg;
464: /* PetscOptionsBegin/End are called in KSPSetFromOptions. They
465: don't need to be called here.
466: */
467: PetscOptionsHeadBegin(PetscOptionsObject, "KSP BiCGStab(L) Options");
469: /* Set number of search directions */
470: PetscOptionsInt("-ksp_bcgsl_ell", "Number of Krylov search directions", "KSPBCGSLSetEll", bcgsl->ell, &this_ell, &flg);
471: if (flg) KSPBCGSLSetEll(ksp, this_ell);
473: /* Set polynomial type */
474: PetscOptionsBool("-ksp_bcgsl_cxpoly", "Polynomial part of BiCGStabL is MinRes + OR", "KSPBCGSLSetPol", flga, &flga, NULL);
475: if (flga) {
476: KSPBCGSLSetPol(ksp, PETSC_TRUE);
477: } else {
478: flg = PETSC_FALSE;
479: PetscOptionsBool("-ksp_bcgsl_mrpoly", "Polynomial part of BiCGStabL is MinRes", "KSPBCGSLSetPol", flg, &flg, NULL);
480: KSPBCGSLSetPol(ksp, PETSC_FALSE);
481: }
483: /* Will computed residual be refreshed? */
484: PetscOptionsReal("-ksp_bcgsl_xres", "Threshold used to decide when to refresh computed residuals", "KSPBCGSLSetXRes", bcgsl->delta, &delta, &flg);
485: if (flg) KSPBCGSLSetXRes(ksp, delta);
487: /* Use pseudoinverse? */
488: flg = bcgsl->pinv;
489: PetscOptionsBool("-ksp_bcgsl_pinv", "Polynomial correction via pseudoinverse", "KSPBCGSLSetUsePseudoinverse", flg, &flg, NULL);
490: KSPBCGSLSetUsePseudoinverse(ksp, flg);
491: PetscOptionsHeadEnd();
492: return 0;
493: }
495: PetscErrorCode KSPSetUp_BCGSL(KSP ksp)
496: {
497: KSP_BCGSL *bcgsl = (KSP_BCGSL *)ksp->data;
498: PetscInt ell = bcgsl->ell, ldMZ = ell + 1;
500: KSPSetWorkVecs(ksp, 6 + 2 * ell);
501: PetscMalloc5(ldMZ, &AY0c, ldMZ, &AYlc, ldMZ, &AYtc, ldMZ * ldMZ, &MZa, ldMZ * ldMZ, &MZb);
502: PetscBLASIntCast(5 * ell, &bcgsl->lwork);
503: PetscMalloc5(bcgsl->lwork, &bcgsl->work, ell, &bcgsl->s, ell * ell, &bcgsl->u, ell * ell, &bcgsl->v, 5 * ell, &bcgsl->realwork);
504: return 0;
505: }
507: PetscErrorCode KSPReset_BCGSL(KSP ksp)
508: {
509: KSP_BCGSL *bcgsl = (KSP_BCGSL *)ksp->data;
511: VecDestroyVecs(ksp->nwork, &ksp->work);
512: PetscFree5(AY0c, AYlc, AYtc, MZa, MZb);
513: PetscFree5(bcgsl->work, bcgsl->s, bcgsl->u, bcgsl->v, bcgsl->realwork);
514: return 0;
515: }
517: PetscErrorCode KSPDestroy_BCGSL(KSP ksp)
518: {
519: KSPReset_BCGSL(ksp);
520: KSPDestroyDefault(ksp);
521: return 0;
522: }
524: /*MC
525: KSPBCGSL - Implements a slight variant of the Enhanced
526: BiCGStab(L) algorithm in (3) and (2). The variation
527: concerns cases when either kappa0**2 or kappa1**2 is
528: negative due to round-off. Kappa0 has also been pulled
529: out of the denominator in the formula for ghat.
531: References:
532: + * - G.L.G. Sleijpen, H.A. van der Vorst, "An overview of
533: approaches for the stable computation of hybrid BiCG
534: methods", Applied Numerical Mathematics: Transactions
535: f IMACS, 19(3), 1996.
536: . * - G.L.G. Sleijpen, H.A. van der Vorst, D.R. Fokkema,
537: "BiCGStab(L) and other hybrid BiCG methods",
538: Numerical Algorithms, 7, 1994.
539: - * - D.R. Fokkema, "Enhanced implementation of BiCGStab(L)
540: for solving linear systems of equations", preprint
541: from www.citeseer.com.
543: Contributed by: Joel M. Malard, email jm.malard@pnl.gov
545: Options Database Keys:
546: + -ksp_bcgsl_ell <ell> Number of Krylov search directions, defaults to 2 -- KSPBCGSLSetEll()
547: . -ksp_bcgsl_cxpol - Use a convex function of the MinRes and OR polynomials after the BiCG step instead of default MinRes -- KSPBCGSLSetPol()
548: . -ksp_bcgsl_mrpoly - Use the default MinRes polynomial after the BiCG step -- KSPBCGSLSetPol()
549: . -ksp_bcgsl_xres <res> Threshold used to decide when to refresh computed residuals -- KSPBCGSLSetXRes()
550: - -ksp_bcgsl_pinv <true/false> - (de)activate use of pseudoinverse -- KSPBCGSLSetUsePseudoinverse()
552: Level: beginner
554: .seealso: `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPFGMRES`, `KSPBCGS`, `KSPSetPCSide()`, `KSPBCGSLSetEll()`, `KSPBCGSLSetXRes()`
556: M*/
557: PETSC_EXTERN PetscErrorCode KSPCreate_BCGSL(KSP ksp)
558: {
559: KSP_BCGSL *bcgsl;
561: /* allocate BiCGStab(L) context */
562: PetscNew(&bcgsl);
563: ksp->data = (void *)bcgsl;
565: KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 3);
566: KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_RIGHT, 2);
567: KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_RIGHT, 1);
569: ksp->ops->setup = KSPSetUp_BCGSL;
570: ksp->ops->solve = KSPSolve_BCGSL;
571: ksp->ops->reset = KSPReset_BCGSL;
572: ksp->ops->destroy = KSPDestroy_BCGSL;
573: ksp->ops->buildsolution = KSPBuildSolutionDefault;
574: ksp->ops->buildresidual = KSPBuildResidualDefault;
575: ksp->ops->setfromoptions = KSPSetFromOptions_BCGSL;
576: ksp->ops->view = KSPView_BCGSL;
578: /* Let the user redefine the number of directions vectors */
579: bcgsl->ell = 2;
581: /*Choose between a single MR step or an averaged MR/OR */
582: bcgsl->bConvex = PETSC_FALSE;
584: bcgsl->pinv = PETSC_TRUE; /* Use the reliable method by default */
586: /* Set the threshold for when exact residuals will be used */
587: bcgsl->delta = 0.0;
588: return 0;
589: }