Actual source code: minres.c
1: #include "src/ksp/ksp/kspimpl.h"
3: typedef struct {
4: PetscReal haptol;
5: } KSP_MINRES;
9: PetscErrorCode KSPSetUp_MINRES(KSP ksp)
10: {
14: if (ksp->pc_side == PC_RIGHT) {
15: SETERRQ(2,"No right preconditioning for KSPMINRES");
16: } else if (ksp->pc_side == PC_SYMMETRIC) {
17: SETERRQ(2,"No symmetric preconditioning for KSPMINRES");
18: }
19: KSPDefaultGetWork(ksp,9);
20: return(0);
21: }
26: PetscErrorCode KSPSolve_MINRES(KSP ksp)
27: {
29: PetscInt i;
30: PetscScalar alpha,malpha,beta,mbeta,ibeta,betaold,eta,c=1.0,ceta,cold=1.0,coold,s=0.0,sold=0.0,soold;
31: PetscScalar rho0,rho1,irho1,rho2,mrho2,rho3,mrho3,mone = -1.0,zero = 0.0,dp = 0.0;
32: PetscReal np;
33: Vec X,B,R,Z,U,V,W,UOLD,VOLD,WOLD,WOOLD;
34: Mat Amat,Pmat;
35: MatStructure pflag;
36: KSP_MINRES *minres = (KSP_MINRES*)ksp->data;
37: PetscTruth diagonalscale;
40: PCDiagonalScale(ksp->pc,&diagonalscale);
41: if (diagonalscale) SETERRQ1(PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",ksp->type_name);
43: X = ksp->vec_sol;
44: B = ksp->vec_rhs;
45: R = ksp->work[0];
46: Z = ksp->work[1];
47: U = ksp->work[2];
48: V = ksp->work[3];
49: W = ksp->work[4];
50: UOLD = ksp->work[5];
51: VOLD = ksp->work[6];
52: WOLD = ksp->work[7];
53: WOOLD = ksp->work[8];
55: PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);
57: ksp->its = 0;
59: VecSet(&zero,UOLD); /* u_old <- 0 */
60: VecCopy(UOLD,VOLD); /* v_old <- 0 */
61: VecCopy(UOLD,W); /* w <- 0 */
62: VecCopy(UOLD,WOLD); /* w_old <- 0 */
64: if (!ksp->guess_zero) {
65: KSP_MatMult(ksp,Amat,X,R); /* r <- b - A*x */
66: VecAYPX(&mone,B,R);
67: } else {
68: VecCopy(B,R); /* r <- b (x is 0) */
69: }
71: KSP_PCApply(ksp,R,Z); /* z <- B*r */
73: VecDot(R,Z,&dp);
74: if (PetscAbsScalar(dp) < minres->haptol) {
75: PetscLogInfo(ksp,"KSPSolve_MINRES:Detected happy breakdown %g tolerance %g\n",PetscAbsScalar(dp),minres->haptol);
76: dp = PetscAbsScalar(dp); /* tiny number, can't use 0.0, cause divided by below */
77: if (dp == 0.0) {
78: ksp->reason = KSP_CONVERGED_ATOL;
79: return(0);
80: }
81: }
83: #if !defined(PETSC_USE_COMPLEX)
84: if (dp < 0.0) {
85: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
86: return(0);
87: }
88: #endif
89: dp = PetscSqrtScalar(dp);
90: beta = dp; /* beta <- sqrt(r'*z */
91: eta = beta;
93: VecCopy(R,V);
94: VecCopy(Z,U);
95: ibeta = 1.0 / beta;
96: VecScale(&ibeta,V); /* v <- r / beta */
97: VecScale(&ibeta,U); /* u <- z / beta */
99: VecNorm(Z,NORM_2,&np); /* np <- ||z|| */
101: KSPLogResidualHistory(ksp,np);
102: KSPMonitor(ksp,0,np); /* call any registered monitor routines */
103: ksp->rnorm = np;
104: (*ksp->converged)(ksp,0,np,&ksp->reason,ksp->cnvP); /* test for convergence */
105: if (ksp->reason) return(0);
107: i = 0;
108: do {
109: ksp->its = i+1;
111: /* Lanczos */
113: KSP_MatMult(ksp,Amat,U,R); /* r <- A*u */
114: VecDot(U,R,&alpha); /* alpha <- r'*u */
115: KSP_PCApply(ksp,R,Z); /* z <- B*r */
117: malpha = - alpha;
118: VecAXPY(&malpha,V,R); /* r <- r - alpha v */
119: VecAXPY(&malpha,U,Z); /* z <- z - alpha u */
120: mbeta = - beta;
121: VecAXPY(&mbeta,VOLD,R); /* r <- r - beta v_old */
122: VecAXPY(&mbeta,UOLD,Z); /* z <- z - beta u_old */
124: betaold = beta;
126: VecDot(R,Z,&dp);
127: if (PetscAbsScalar(dp) < minres->haptol) {
128: PetscLogInfo(ksp,"KSPSolve_MINRES:Detected happy breakdown %g tolerance %g\n",PetscAbsScalar(dp),minres->haptol);
129: dp = PetscAbsScalar(dp); /* tiny number, can we use 0.0? */
130: }
132: #if !defined(PETSC_USE_COMPLEX)
133: if (dp < 0.0) {
134: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
135: break;
136: }
138: #endif
139: beta = PetscSqrtScalar(dp); /* beta <- sqrt(r'*z) */
141: /* QR factorisation */
143: coold = cold; cold = c; soold = sold; sold = s;
145: rho0 = cold * alpha - coold * sold * betaold;
146: rho1 = PetscSqrtScalar(rho0*rho0 + beta*beta);
147: rho2 = sold * alpha + coold * cold * betaold;
148: rho3 = soold * betaold;
150: /* Givens rotation */
152: c = rho0 / rho1;
153: s = beta / rho1;
155: /* Update */
157: VecCopy(WOLD,WOOLD); /* w_oold <- w_old */
158: VecCopy(W,WOLD); /* w_old <- w */
159:
160: VecCopy(U,W); /* w <- u */
161: mrho2 = - rho2;
162: VecAXPY(&mrho2,WOLD,W); /* w <- w - rho2 w_old */
163: mrho3 = - rho3;
164: VecAXPY(&mrho3,WOOLD,W); /* w <- w - rho3 w_oold */
165: irho1 = 1.0 / rho1;
166: VecScale(&irho1,W); /* w <- w / rho1 */
168: ceta = c * eta;
169: VecAXPY(&ceta,W,X); /* x <- x + c eta w */
170: eta = - s * eta;
172: VecCopy(V,VOLD);
173: VecCopy(U,UOLD);
174: VecCopy(R,V);
175: VecCopy(Z,U);
176: ibeta = 1.0 / beta;
177: VecScale(&ibeta,V); /* v <- r / beta */
178: VecScale(&ibeta,U); /* u <- z / beta */
179:
180: np = ksp->rnorm * PetscAbsScalar(s);
182: ksp->rnorm = np;
183: KSPLogResidualHistory(ksp,np);
184: KSPMonitor(ksp,i+1,np);
185: (*ksp->converged)(ksp,i+1,np,&ksp->reason,ksp->cnvP); /* test for convergence */
186: if (ksp->reason) break;
187: i++;
188: } while (i<ksp->max_it);
189: if (i >= ksp->max_it) {
190: ksp->reason = KSP_DIVERGED_ITS;
191: }
192: return(0);
193: }
195: /*MC
196: KSPMINRES - This code implements the MINRES (Minimum Residual) method.
198: Options Database Keys:
199: . see KSPSolve()
201: Level: beginner
203: Contributed by: Robert Scheichl: maprs@maths.bath.ac.uk
205: Notes: The operator and the preconditioner must be positive definite for this method
206: Reference: Paige & Saunders, 1975.
208: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPCG
209: M*/
213: PetscErrorCode KSPCreate_MINRES(KSP ksp)
214: {
215: KSP_MINRES *minres;
220: ksp->pc_side = PC_LEFT;
221: PetscNew(KSP_MINRES,&minres);
222: minres->haptol = 1.e-18;
223: ksp->data = (void*)minres;
225: /*
226: Sets the functions that are associated with this data structure
227: (in C++ this is the same as defining virtual functions)
228: */
229: ksp->ops->setup = KSPSetUp_MINRES;
230: ksp->ops->solve = KSPSolve_MINRES;
231: ksp->ops->destroy = KSPDefaultDestroy;
232: ksp->ops->setfromoptions = 0;
233: ksp->ops->buildsolution = KSPDefaultBuildSolution;
234: ksp->ops->buildresidual = KSPDefaultBuildResidual;
235: return(0);
236: }