Actual source code: minres.c
2: #include <petsc/private/kspimpl.h>
4: typedef struct {
5: PetscReal haptol;
6: } KSP_MINRES;
8: static PetscErrorCode KSPSetUp_MINRES(KSP ksp)
9: {
10: PetscFunctionBegin;
11: PetscCall(KSPSetWorkVecs(ksp, 9));
12: PetscFunctionReturn(PETSC_SUCCESS);
13: }
15: #define KSPMinresSwap3(V1, V2, V3) \
16: do { \
17: Vec T = V1; \
18: V1 = V2; \
19: V2 = V3; \
20: V3 = T; \
21: } while (0)
23: static PetscErrorCode KSPSolve_MINRES(KSP ksp)
24: {
25: PetscInt i;
26: PetscScalar alpha, beta, betaold, eta, c = 1.0, ceta, cold = 1.0, coold, s = 0.0, sold = 0.0, soold;
27: PetscScalar rho0, rho1, rho2, rho3, dp = 0.0;
28: const PetscScalar none = -1.0;
29: PetscReal np;
30: Vec X, B, R, Z, U, V, W, UOLD, VOLD, WOLD, WOOLD;
31: Mat Amat;
32: KSP_MINRES *minres = (KSP_MINRES *)ksp->data;
33: PetscBool diagonalscale;
35: PetscFunctionBegin;
36: PetscCall(PCGetDiagonalScale(ksp->pc, &diagonalscale));
37: PetscCheck(!diagonalscale, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Krylov method %s does not support diagonal scaling", ((PetscObject)ksp)->type_name);
39: X = ksp->vec_sol;
40: B = ksp->vec_rhs;
41: R = ksp->work[0];
42: Z = ksp->work[1];
43: U = ksp->work[2];
44: V = ksp->work[3];
45: W = ksp->work[4];
46: UOLD = ksp->work[5];
47: VOLD = ksp->work[6];
48: WOLD = ksp->work[7];
49: WOOLD = ksp->work[8];
51: PetscCall(PCGetOperators(ksp->pc, &Amat, NULL));
53: ksp->its = 0;
55: if (!ksp->guess_zero) {
56: PetscCall(KSP_MatMult(ksp, Amat, X, R)); /* r <- b - A*x */
57: PetscCall(VecAYPX(R, -1.0, B));
58: } else {
59: PetscCall(VecCopy(B, R)); /* r <- b (x is 0) */
60: }
61: PetscCall(KSP_PCApply(ksp, R, Z)); /* z <- B*r */
62: PetscCall(VecNorm(Z, NORM_2, &np)); /* np <- ||z|| */
63: KSPCheckNorm(ksp, np);
64: PetscCall(VecDot(R, Z, &dp));
65: KSPCheckDot(ksp, dp);
67: if (PetscRealPart(dp) < minres->haptol && np > minres->haptol) {
68: PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_CONV_FAILED, "Detected indefinite operator %g tolerance %g", (double)PetscRealPart(dp), (double)minres->haptol);
69: PetscCall(PetscInfo(ksp, "Detected indefinite operator %g tolerance %g\n", (double)PetscRealPart(dp), (double)minres->haptol));
70: ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
71: PetscFunctionReturn(PETSC_SUCCESS);
72: }
74: ksp->rnorm = 0.0;
75: if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = np;
76: PetscCall(KSPLogResidualHistory(ksp, ksp->rnorm));
77: PetscCall(KSPMonitor(ksp, 0, ksp->rnorm));
78: PetscCall((*ksp->converged)(ksp, 0, ksp->rnorm, &ksp->reason, ksp->cnvP)); /* test for convergence */
79: if (ksp->reason) PetscFunctionReturn(PETSC_SUCCESS);
81: dp = PetscAbsScalar(dp);
82: dp = PetscSqrtScalar(dp);
83: beta = dp; /* beta <- sqrt(r'*z */
84: eta = beta;
85: PetscCall(VecAXPBY(V, 1.0 / beta, 0, R)); /* v <- r / beta */
86: PetscCall(VecAXPBY(U, 1.0 / beta, 0, Z)); /* u <- z / beta */
88: i = 0;
89: do {
90: ksp->its = i + 1;
92: /* Lanczos */
94: PetscCall(KSP_MatMult(ksp, Amat, U, R)); /* r <- A*u */
95: PetscCall(VecDot(U, R, &alpha)); /* alpha <- r'*u */
96: PetscCall(KSP_PCApply(ksp, R, Z)); /* z <- B*r */
98: if (ksp->its > 1) {
99: Vec T[2];
100: PetscScalar alphas[] = {-alpha, -beta};
101: /* r <- r - alpha v - beta v_old */
102: T[0] = V;
103: T[1] = VOLD;
104: PetscCall(VecMAXPY(R, 2, alphas, T));
105: /* z <- z - alpha u - beta u_old */
106: T[0] = U;
107: T[1] = UOLD;
108: PetscCall(VecMAXPY(Z, 2, alphas, T));
109: } else {
110: PetscCall(VecAXPY(R, -alpha, V)); /* r <- r - alpha v */
111: PetscCall(VecAXPY(Z, -alpha, U)); /* z <- z - alpha u */
112: }
114: betaold = beta;
116: PetscCall(VecDot(R, Z, &dp));
117: KSPCheckDot(ksp, dp);
118: dp = PetscAbsScalar(dp);
119: beta = PetscSqrtScalar(dp); /* beta <- sqrt(r'*z) */
121: /* QR factorisation */
123: coold = cold;
124: cold = c;
125: soold = sold;
126: sold = s;
128: rho0 = cold * alpha - coold * sold * betaold;
129: rho1 = PetscSqrtScalar(rho0 * rho0 + beta * beta);
130: rho2 = sold * alpha + coold * cold * betaold;
131: rho3 = soold * betaold;
133: /* Givens rotation */
135: c = rho0 / rho1;
136: s = beta / rho1;
138: /* Update */
139: /* w_oold <- w_old */
140: /* w_old <- w */
141: KSPMinresSwap3(WOOLD, WOLD, W);
143: /* w <- (u - rho2 w_old - rho3 w_oold)/rho1 */
144: PetscCall(VecAXPBY(W, 1.0 / rho1, 0.0, U));
145: if (ksp->its > 1) {
146: Vec T[] = {WOLD, WOOLD};
147: PetscScalar alphas[] = {-rho2 / rho1, -rho3 / rho1};
148: PetscInt nv = (ksp->its == 2 ? 1 : 2);
150: PetscCall(VecMAXPY(W, nv, alphas, T));
151: }
153: ceta = c * eta;
154: PetscCall(VecAXPY(X, ceta, W)); /* x <- x + c eta w */
156: /*
157: when dp is really small we have either convergence or an indefinite operator so compute true
158: residual norm to check for convergence
159: */
160: if (PetscRealPart(dp) < minres->haptol) {
161: PetscCall(PetscInfo(ksp, "Possible indefinite operator %g tolerance %g\n", (double)PetscRealPart(dp), (double)minres->haptol));
162: PetscCall(KSP_MatMult(ksp, Amat, X, VOLD));
163: PetscCall(VecAXPY(VOLD, none, B));
164: PetscCall(VecNorm(VOLD, NORM_2, &np));
165: KSPCheckNorm(ksp, np);
166: } else {
167: /* otherwise compute new residual norm via recurrence relation */
168: np *= PetscAbsScalar(s);
169: }
171: if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = np;
172: PetscCall(KSPLogResidualHistory(ksp, ksp->rnorm));
173: PetscCall(KSPMonitor(ksp, i + 1, ksp->rnorm));
174: PetscCall((*ksp->converged)(ksp, i + 1, ksp->rnorm, &ksp->reason, ksp->cnvP)); /* test for convergence */
175: if (ksp->reason) break;
177: if (PetscRealPart(dp) < minres->haptol) {
178: PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_CONV_FAILED, "Detected indefinite operator %g tolerance %g", (double)PetscRealPart(dp), (double)minres->haptol);
179: PetscCall(PetscInfo(ksp, "Detected indefinite operator %g tolerance %g\n", (double)PetscRealPart(dp), (double)minres->haptol));
180: ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
181: break;
182: }
184: eta = -s * eta;
185: KSPMinresSwap3(VOLD, V, R);
186: KSPMinresSwap3(UOLD, U, Z);
187: PetscCall(VecScale(V, 1.0 / beta)); /* v <- r / beta */
188: PetscCall(VecScale(U, 1.0 / beta)); /* u <- z / beta */
190: i++;
191: } while (i < ksp->max_it);
192: if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
193: PetscFunctionReturn(PETSC_SUCCESS);
194: }
196: /*MC
197: KSPMINRES - This code implements the MINRES (Minimum Residual) method.
199: Level: beginner
201: Notes:
202: The operator and the preconditioner must be symmetric and the preconditioner must be positive definite for this method.
204: Supports only left preconditioning.
206: Reference:
207: . * - Paige & Saunders, Solution of sparse indefinite systems of linear equations, SIAM J. Numer. Anal. 12, 1975.
209: Contributed by:
210: Robert Scheichl: maprs@maths.bath.ac.uk
212: .seealso: [](chapter_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPCG`, `KSPCR`
213: M*/
214: PETSC_EXTERN PetscErrorCode KSPCreate_MINRES(KSP ksp)
215: {
216: KSP_MINRES *minres;
218: PetscFunctionBegin;
219: PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 3));
220: PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_LEFT, 1));
221: PetscCall(PetscNew(&minres));
223: /* this parameter is arbitrary; but e-50 didn't work for __float128 in one example */
224: #if defined(PETSC_USE_REAL___FLOAT128)
225: minres->haptol = 1.e-100;
226: #elif defined(PETSC_USE_REAL_SINGLE)
227: minres->haptol = 1.e-25;
228: #else
229: minres->haptol = 1.e-50;
230: #endif
231: ksp->data = (void *)minres;
233: /*
234: Sets the functions that are associated with this data structure
235: (in C++ this is the same as defining virtual functions)
236: */
237: ksp->ops->setup = KSPSetUp_MINRES;
238: ksp->ops->solve = KSPSolve_MINRES;
239: ksp->ops->destroy = KSPDestroyDefault;
240: ksp->ops->setfromoptions = NULL;
241: ksp->ops->buildsolution = KSPBuildSolutionDefault;
242: ksp->ops->buildresidual = KSPBuildResidualDefault;
243: PetscFunctionReturn(PETSC_SUCCESS);
244: }