Actual source code: minres.c
petsc-3.7.3 2016-08-01
2: #include <petsc/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(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"No right preconditioning for KSPMINRES");
16: else if (ksp->pc_side == PC_SYMMETRIC) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"No symmetric preconditioning for KSPMINRES");
17: KSPSetWorkVecs(ksp,9);
18: return(0);
19: }
24: PetscErrorCode KSPSolve_MINRES(KSP ksp)
25: {
26: PetscErrorCode ierr;
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,rho3,dp = 0.0;
30: const PetscScalar none = -1.0;
31: PetscReal np;
32: Vec X,B,R,Z,U,V,W,UOLD,VOLD,WOLD,WOOLD;
33: Mat Amat,Pmat;
34: KSP_MINRES *minres = (KSP_MINRES*)ksp->data;
35: PetscBool diagonalscale;
38: PCGetDiagonalScale(ksp->pc,&diagonalscale);
39: if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),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);
55: ksp->its = 0;
57: VecSet(UOLD,0.0); /* u_old <- 0 */
58: VecSet(VOLD,0.0); /* v_old <- 0 */
59: VecSet(W,0.0); /* w <- 0 */
60: VecSet(WOLD,0.0); /* 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: }
68: KSP_PCApply(ksp,R,Z); /* z <- B*r */
69: VecNorm(Z,NORM_2,&np); /* np <- ||z|| */
70: VecDot(R,Z,&dp);
72: if (PetscRealPart(dp) < minres->haptol && np > minres->haptol) {
73: if (ksp->errorifnotconverged) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_CONV_FAILED,"Detected indefinite operator %g tolerance %g",(double)PetscRealPart(dp),(double)minres->haptol);
74: PetscInfo2(ksp,"Detected indefinite operator %g tolerance %g\n",(double)PetscRealPart(dp),(double)minres->haptol);
75: ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
76: return(0);
77: }
79: KSPLogResidualHistory(ksp,np);
80: KSPMonitor(ksp,0,np);
81: ksp->rnorm = np;
82: (*ksp->converged)(ksp,0,np,&ksp->reason,ksp->cnvP); /* test for convergence */
83: if (ksp->reason) return(0);
85: dp = PetscAbsScalar(dp);
86: dp = PetscSqrtScalar(dp);
87: beta = dp; /* beta <- sqrt(r'*z */
88: eta = beta;
90: VecCopy(R,V);
91: VecCopy(Z,U);
92: ibeta = 1.0 / beta;
93: VecScale(V,ibeta); /* v <- r / beta */
94: VecScale(U,ibeta); /* u <- z / beta */
96: i = 0;
97: do {
98: ksp->its = i+1;
100: /* Lanczos */
102: KSP_MatMult(ksp,Amat,U,R); /* r <- A*u */
103: VecDot(U,R,&alpha); /* alpha <- r'*u */
104: KSP_PCApply(ksp,R,Z); /* z <- B*r */
106: VecAXPY(R,-alpha,V); /* r <- r - alpha v */
107: VecAXPY(Z,-alpha,U); /* z <- z - alpha u */
108: VecAXPY(R,-beta,VOLD); /* r <- r - beta v_old */
109: VecAXPY(Z,-beta,UOLD); /* z <- z - beta u_old */
111: betaold = beta;
113: VecDot(R,Z,&dp);
114: dp = PetscAbsScalar(dp);
115: beta = PetscSqrtScalar(dp); /* beta <- sqrt(r'*z) */
117: /* QR factorisation */
119: coold = cold; cold = c; soold = sold; sold = s;
121: rho0 = cold * alpha - coold * sold * betaold;
122: rho1 = PetscSqrtScalar(rho0*rho0 + beta*beta);
123: rho2 = sold * alpha + coold * cold * betaold;
124: rho3 = soold * betaold;
126: /* Givens rotation */
128: c = rho0 / rho1;
129: s = beta / rho1;
131: /* Update */
133: VecCopy(WOLD,WOOLD); /* w_oold <- w_old */
134: VecCopy(W,WOLD); /* w_old <- w */
136: VecCopy(U,W); /* w <- u */
137: VecAXPY(W,-rho2,WOLD); /* w <- w - rho2 w_old */
138: VecAXPY(W,-rho3,WOOLD); /* w <- w - rho3 w_oold */
139: irho1 = 1.0 / rho1;
140: VecScale(W,irho1); /* w <- w / rho1 */
142: ceta = c * eta;
143: VecAXPY(X,ceta,W); /* x <- x + c eta w */
145: /*
146: when dp is really small we have either convergence or an indefinite operator so compute true
147: residual norm to check for convergence
148: */
149: if (PetscRealPart(dp) < minres->haptol) {
150: PetscInfo2(ksp,"Possible indefinite operator %g tolerance %g\n",(double)PetscRealPart(dp),(double)minres->haptol);
151: KSP_MatMult(ksp,Amat,X,VOLD);
152: VecAXPY(VOLD,none,B);
153: VecNorm(VOLD,NORM_2,&np);
154: } else {
155: /* otherwise compute new residual norm via recurrence relation */
156: np = ksp->rnorm * PetscAbsScalar(s);
157: }
159: ksp->rnorm = np;
160: KSPLogResidualHistory(ksp,np);
161: KSPMonitor(ksp,i+1,np);
162: (*ksp->converged)(ksp,i+1,np,&ksp->reason,ksp->cnvP); /* test for convergence */
163: if (ksp->reason) break;
165: if (PetscRealPart(dp) < minres->haptol) {
166: if (ksp->errorifnotconverged) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_CONV_FAILED,"Detected indefinite operator %g tolerance %g",(double)PetscRealPart(dp),(double)minres->haptol);
167: PetscInfo2(ksp,"Detected indefinite operator %g tolerance %g\n",(double)PetscRealPart(dp),(double)minres->haptol);
168: ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
169: break;
170: }
172: eta = -s * eta;
173: VecCopy(V,VOLD);
174: VecCopy(U,UOLD);
175: VecCopy(R,V);
176: VecCopy(Z,U);
177: ibeta = 1.0 / beta;
178: VecScale(V,ibeta); /* v <- r / beta */
179: VecScale(U,ibeta); /* u <- z / beta */
181: i++;
182: } while (i<ksp->max_it);
183: if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
184: return(0);
185: }
187: /*MC
188: KSPMINRES - This code implements the MINRES (Minimum Residual) method.
190: Options Database Keys:
191: . see KSPSolve()
193: Level: beginner
195: Notes: The operator and the preconditioner must be symmetric and the preconditioner must
196: be positive definite for this method.
197: Supports only left preconditioning.
199: Reference: Paige & Saunders, 1975.
201: Contributed by: Robert Scheichl: maprs@maths.bath.ac.uk
203: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPCG, KSPCR
204: M*/
207: PETSC_EXTERN PetscErrorCode KSPCreate_MINRES(KSP ksp)
208: {
209: KSP_MINRES *minres;
213: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
214: PetscNewLog(ksp,&minres);
216: /* this parameter is arbitrary; but e-50 didn't work for __float128 in one example */
217: #if defined(PETSC_USE_REAL___FLOAT128)
218: minres->haptol = 1.e-100;
219: #else
220: minres->haptol = 1.e-50;
221: #endif
222: ksp->data = (void*)minres;
224: /*
225: Sets the functions that are associated with this data structure
226: (in C++ this is the same as defining virtual functions)
227: */
228: ksp->ops->setup = KSPSetUp_MINRES;
229: ksp->ops->solve = KSPSolve_MINRES;
230: ksp->ops->destroy = KSPDestroyDefault;
231: ksp->ops->setfromoptions = 0;
232: ksp->ops->buildsolution = KSPBuildSolutionDefault;
233: ksp->ops->buildresidual = KSPBuildResidualDefault;
234: return(0);
235: }