Actual source code: minres.c
2: #include <private/kspimpl.h>
4: typedef struct {
5: PetscReal haptol;
6: } KSP_MINRES;
10: PetscErrorCode KSPSetUp_MINRES(KSP ksp)
11: {
15: if (ksp->pc_side == PC_RIGHT) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"No right preconditioning for KSPMINRES");
16: else if (ksp->pc_side == PC_SYMMETRIC) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"No symmetric preconditioning for KSPMINRES");
17: KSPDefaultGetWork(ksp,9);
18: return(0);
19: }
24: PetscErrorCode KSPSolve_MINRES(KSP ksp)
25: {
27: PetscInt i;
28: PetscScalar alpha,beta,ibeta,betaold,eta,c=1.0,ceta,cold=1.0,coold,s=0.0,sold=0.0,soold;
29: PetscScalar rho0,rho1,irho1,rho2,mrho2,rho3,mrho3,dp = 0.0;
30: PetscReal np;
31: Vec X,B,R,Z,U,V,W,UOLD,VOLD,WOLD,WOOLD;
32: Mat Amat,Pmat;
33: MatStructure pflag;
34: KSP_MINRES *minres = (KSP_MINRES*)ksp->data;
35: PetscBool diagonalscale;
38: PCGetDiagonalScale(ksp->pc,&diagonalscale);
39: if (diagonalscale) SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
41: X = ksp->vec_sol;
42: B = ksp->vec_rhs;
43: R = ksp->work[0];
44: Z = ksp->work[1];
45: U = ksp->work[2];
46: V = ksp->work[3];
47: W = ksp->work[4];
48: UOLD = ksp->work[5];
49: VOLD = ksp->work[6];
50: WOLD = ksp->work[7];
51: WOOLD = ksp->work[8];
53: PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);
55: ksp->its = 0;
57: VecSet(UOLD,0.0); /* u_old <- 0 */
58: VecCopy(UOLD,VOLD); /* v_old <- 0 */
59: VecCopy(UOLD,W); /* w <- 0 */
60: VecCopy(UOLD,WOLD); /* w_old <- 0 */
62: if (!ksp->guess_zero) {
63: KSP_MatMult(ksp,Amat,X,R); /* r <- b - A*x */
64: VecAYPX(R,-1.0,B);
65: } else {
66: VecCopy(B,R); /* r <- b (x is 0) */
67: }
69: KSP_PCApply(ksp,R,Z); /* z <- B*r */
71: VecDot(R,Z,&dp);
72: if (PetscAbsScalar(dp) < minres->haptol) {
73: PetscInfo2(ksp,"Detected happy breakdown %G tolerance %G\n",PetscAbsScalar(dp),minres->haptol);
74: dp = PetscAbsScalar(dp); /* tiny number, can't use 0.0, cause divided by below */
75: if (dp == 0.0) {
76: ksp->reason = KSP_CONVERGED_ATOL;
77: return(0);
78: }
79: }
81: #if !defined(PETSC_USE_COMPLEX)
82: if (dp < 0.0) {
83: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
84: return(0);
85: }
86: #endif
87: dp = PetscSqrtScalar(dp);
88: beta = dp; /* beta <- sqrt(r'*z */
89: eta = beta;
91: VecCopy(R,V);
92: VecCopy(Z,U);
93: ibeta = 1.0 / beta;
94: VecScale(V,ibeta); /* v <- r / beta */
95: VecScale(U,ibeta); /* u <- z / beta */
97: VecNorm(Z,NORM_2,&np); /* np <- ||z|| */
99: KSPLogResidualHistory(ksp,np);
100: KSPMonitor(ksp,0,np);
101: ksp->rnorm = np;
102: (*ksp->converged)(ksp,0,np,&ksp->reason,ksp->cnvP); /* test for convergence */
103: if (ksp->reason) return(0);
105: i = 0;
106: do {
107: ksp->its = i+1;
109: /* Lanczos */
111: KSP_MatMult(ksp,Amat,U,R); /* r <- A*u */
112: VecDot(U,R,&alpha); /* alpha <- r'*u */
113: KSP_PCApply(ksp,R,Z); /* z <- B*r */
115: VecAXPY(R,-alpha,V); /* r <- r - alpha v */
116: VecAXPY(Z,-alpha,U); /* z <- z - alpha u */
117: VecAXPY(R,-beta,VOLD); /* r <- r - beta v_old */
118: VecAXPY(Z,-beta,UOLD); /* z <- z - beta u_old */
120: betaold = beta;
122: VecDot(R,Z,&dp);
123: if (PetscAbsScalar(dp) < minres->haptol) {
124: PetscInfo2(ksp,"Detected happy breakdown %G tolerance %G\n",PetscAbsScalar(dp),minres->haptol);
125: dp = PetscAbsScalar(dp); /* tiny number, can we use 0.0? */
126: }
128: #if !defined(PETSC_USE_COMPLEX)
129: if (dp < 0.0) {
130: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
131: break;
132: }
134: #endif
135: beta = PetscSqrtScalar(dp); /* beta <- sqrt(r'*z) */
137: /* QR factorisation */
139: coold = cold; cold = c; soold = sold; sold = s;
141: rho0 = cold * alpha - coold * sold * betaold;
142: rho1 = PetscSqrtScalar(rho0*rho0 + beta*beta);
143: rho2 = sold * alpha + coold * cold * betaold;
144: rho3 = soold * betaold;
146: /* Givens rotation */
148: c = rho0 / rho1;
149: s = beta / rho1;
151: /* Update */
153: VecCopy(WOLD,WOOLD); /* w_oold <- w_old */
154: VecCopy(W,WOLD); /* w_old <- w */
155:
156: VecCopy(U,W); /* w <- u */
157: mrho2 = - rho2;
158: VecAXPY(W,mrho2,WOLD); /* w <- w - rho2 w_old */
159: mrho3 = - rho3;
160: VecAXPY(W,mrho3,WOOLD); /* w <- w - rho3 w_oold */
161: irho1 = 1.0 / rho1;
162: VecScale(W,irho1); /* w <- w / rho1 */
164: ceta = c * eta;
165: VecAXPY(X,ceta,W); /* x <- x + c eta w */
166: eta = - s * eta;
168: VecCopy(V,VOLD);
169: VecCopy(U,UOLD);
170: VecCopy(R,V);
171: VecCopy(Z,U);
172: ibeta = 1.0 / beta;
173: VecScale(V,ibeta); /* v <- r / beta */
174: VecScale(U,ibeta); /* u <- z / beta */
175:
176: np = ksp->rnorm * PetscAbsScalar(s);
178: ksp->rnorm = np;
179: KSPLogResidualHistory(ksp,np);
180: KSPMonitor(ksp,i+1,np);
181: (*ksp->converged)(ksp,i+1,np,&ksp->reason,ksp->cnvP); /* test for convergence */
182: if (ksp->reason) break;
183: i++;
184: } while (i<ksp->max_it);
185: if (i >= ksp->max_it) {
186: ksp->reason = KSP_DIVERGED_ITS;
187: }
188: return(0);
189: }
191: /*MC
192: KSPMINRES - This code implements the MINRES (Minimum Residual) method.
194: Options Database Keys:
195: . see KSPSolve()
197: Level: beginner
199: Notes: The operator and the preconditioner must be symmetric and the preconditioner must
200: be positive definite for this method.
201: Supports only left preconditioning.
203: Reference: Paige & Saunders, 1975.
205: Contributed by: Robert Scheichl: maprs@maths.bath.ac.uk
207: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPCG, KSPCR
208: M*/
212: PetscErrorCode KSPCreate_MINRES(KSP ksp)
213: {
214: KSP_MINRES *minres;
218: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
219: PetscNewLog(ksp,KSP_MINRES,&minres);
220: minres->haptol = 1.e-18;
221: ksp->data = (void*)minres;
223: /*
224: Sets the functions that are associated with this data structure
225: (in C++ this is the same as defining virtual functions)
226: */
227: ksp->ops->setup = KSPSetUp_MINRES;
228: ksp->ops->solve = KSPSolve_MINRES;
229: ksp->ops->destroy = KSPDefaultDestroy;
230: ksp->ops->setfromoptions = 0;
231: ksp->ops->buildsolution = KSPDefaultBuildSolution;
232: ksp->ops->buildresidual = KSPDefaultBuildResidual;
233: return(0);
234: }