Actual source code: tfqmr.c
1: #include <petsc/private/kspimpl.h>
3: static PetscErrorCode KSPSetUp_TFQMR(KSP ksp)
4: {
5: PetscFunctionBegin;
6: PetscCheck(ksp->pc_side != PC_SYMMETRIC, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "no symmetric preconditioning for KSPTFQMR");
7: PetscCall(KSPSetWorkVecs(ksp, 9));
8: PetscFunctionReturn(PETSC_SUCCESS);
9: }
11: static PetscErrorCode KSPSolve_TFQMR(KSP ksp)
12: {
13: PetscInt i, m;
14: PetscScalar rho, rhoold, a, s, b, eta, etaold, psiold, cf;
15: PetscReal dp, dpold, w, dpest, tau, psi, cm;
16: Vec X, B, V, P, R, RP, T, T1, Q, U, D, AUQ;
18: PetscFunctionBegin;
19: X = ksp->vec_sol;
20: B = ksp->vec_rhs;
21: R = ksp->work[0];
22: RP = ksp->work[1];
23: V = ksp->work[2];
24: T = ksp->work[3];
25: Q = ksp->work[4];
26: P = ksp->work[5];
27: U = ksp->work[6];
28: D = ksp->work[7];
29: T1 = ksp->work[8];
30: AUQ = V;
32: /* Compute initial preconditioned residual */
33: PetscCall(KSPInitialResidual(ksp, X, V, T, R, B));
35: /* Test for nothing to do */
36: PetscCall(VecNorm(R, NORM_2, &dp));
37: KSPCheckNorm(ksp, dp);
38: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
39: if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = dp;
40: else ksp->rnorm = 0.0;
41: ksp->its = 0;
42: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
43: PetscCall(KSPMonitor(ksp, 0, ksp->rnorm));
44: PetscCall((*ksp->converged)(ksp, 0, ksp->rnorm, &ksp->reason, ksp->cnvP));
45: if (ksp->reason) PetscFunctionReturn(PETSC_SUCCESS);
47: /* Make the initial Rp == R */
48: PetscCall(VecCopy(R, RP));
50: /* Set the initial conditions */
51: etaold = 0.0;
52: psiold = 0.0;
53: tau = dp;
54: dpold = dp;
56: PetscCall(VecDot(R, RP, &rhoold)); /* rhoold = (r,rp) */
57: PetscCall(VecCopy(R, U));
58: PetscCall(VecCopy(R, P));
59: PetscCall(KSP_PCApplyBAorAB(ksp, P, V, T));
60: PetscCall(VecSet(D, 0.0));
62: i = 0;
63: do {
64: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
65: ksp->its++;
66: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
67: PetscCall(VecDot(V, RP, &s)); /* s <- (v,rp) */
68: KSPCheckDot(ksp, s);
69: a = rhoold / s; /* a <- rho / s */
70: PetscCall(VecWAXPY(Q, -a, V, U)); /* q <- u - a v */
71: PetscCall(VecWAXPY(T, 1.0, U, Q)); /* t <- u + q */
72: PetscCall(KSP_PCApplyBAorAB(ksp, T, AUQ, T1));
73: PetscCall(VecAXPY(R, -a, AUQ)); /* r <- r - a K (u + q) */
74: PetscCall(VecNorm(R, NORM_2, &dp));
75: KSPCheckNorm(ksp, dp);
76: for (m = 0; m < 2; m++) {
77: if (!m) w = PetscSqrtReal(dp * dpold);
78: else w = dp;
79: psi = w / tau;
80: cm = 1.0 / PetscSqrtReal(1.0 + psi * psi);
81: tau = tau * psi * cm;
82: eta = cm * cm * a;
83: cf = psiold * psiold * etaold / a;
84: if (!m) {
85: PetscCall(VecAYPX(D, cf, U));
86: } else {
87: PetscCall(VecAYPX(D, cf, Q));
88: }
89: PetscCall(VecAXPY(X, eta, D));
91: dpest = PetscSqrtReal(2 * i + m + 2.0) * tau;
92: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
93: if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = dpest;
94: else ksp->rnorm = 0.0;
95: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
96: PetscCall(KSPLogResidualHistory(ksp, ksp->rnorm));
97: PetscCall(KSPMonitor(ksp, i + 1, ksp->rnorm));
98: PetscCall((*ksp->converged)(ksp, i + 1, ksp->rnorm, &ksp->reason, ksp->cnvP));
99: if (ksp->reason) break;
101: etaold = eta;
102: psiold = psi;
103: }
104: if (ksp->reason) break;
106: PetscCall(VecDot(R, RP, &rho)); /* rho <- (r,rp) */
107: b = rho / rhoold; /* b <- rho / rhoold */
108: PetscCall(VecWAXPY(U, b, Q, R)); /* u <- r + b q */
109: PetscCall(VecAXPY(Q, b, P));
110: PetscCall(VecWAXPY(P, b, Q, U)); /* p <- u + b(q + b p) */
111: PetscCall(KSP_PCApplyBAorAB(ksp, P, V, Q)); /* v <- K p */
113: rhoold = rho;
114: dpold = dp;
116: i++;
117: } while (i < ksp->max_it);
118: if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
120: PetscCall(KSPUnwindPreconditioner(ksp, X, T));
121: PetscFunctionReturn(PETSC_SUCCESS);
122: }
124: /*MC
125: KSPTFQMR - A transpose-free QMR (quasi minimal residual) {cite}`f:93`
127: Level: beginner
129: Notes:
130: Supports left and right preconditioning, but not symmetric
132: The "residual norm" computed in this algorithm is actually just an upper bound on the actual residual norm.
133: That is for left preconditioning it is a bound on the preconditioned residual and for right preconditioning
134: it is a bound on the true residual.
136: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPTCQMR`
137: M*/
138: PETSC_EXTERN PetscErrorCode KSPCreate_TFQMR(KSP ksp)
139: {
140: PetscFunctionBegin;
141: PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 3));
142: PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_RIGHT, 2));
143: PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_LEFT, 1));
144: PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_RIGHT, 1));
146: ksp->data = (void *)0;
147: ksp->ops->setup = KSPSetUp_TFQMR;
148: ksp->ops->solve = KSPSolve_TFQMR;
149: ksp->ops->destroy = KSPDestroyDefault;
150: ksp->ops->buildsolution = KSPBuildSolutionDefault;
151: ksp->ops->buildresidual = KSPBuildResidualDefault;
152: ksp->ops->setfromoptions = NULL;
153: ksp->ops->view = NULL;
154: PetscFunctionReturn(PETSC_SUCCESS);
155: }