Actual source code: minres.c
petsc-3.10.5 2019-03-28
2: #include <petsc/private/kspimpl.h>
4: typedef struct {
5: PetscReal haptol;
6: } KSP_MINRES;
8: static PetscErrorCode KSPSetUp_MINRES(KSP ksp)
9: {
13: if (ksp->pc_side == PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"No right preconditioning for KSPMINRES");
14: else if (ksp->pc_side == PC_SYMMETRIC) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"No symmetric preconditioning for KSPMINRES");
15: KSPSetWorkVecs(ksp,9);
16: return(0);
17: }
20: static PetscErrorCode KSPSolve_MINRES(KSP ksp)
21: {
22: PetscErrorCode ierr;
23: PetscInt i;
24: PetscScalar alpha,beta,ibeta,betaold,eta,c=1.0,ceta,cold=1.0,coold,s=0.0,sold=0.0,soold;
25: PetscScalar rho0,rho1,irho1,rho2,rho3,dp = 0.0;
26: const PetscScalar none = -1.0;
27: PetscReal np;
28: Vec X,B,R,Z,U,V,W,UOLD,VOLD,WOLD,WOOLD;
29: Mat Amat,Pmat;
30: KSP_MINRES *minres = (KSP_MINRES*)ksp->data;
31: PetscBool diagonalscale;
34: PCGetDiagonalScale(ksp->pc,&diagonalscale);
35: if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
37: X = ksp->vec_sol;
38: B = ksp->vec_rhs;
39: R = ksp->work[0];
40: Z = ksp->work[1];
41: U = ksp->work[2];
42: V = ksp->work[3];
43: W = ksp->work[4];
44: UOLD = ksp->work[5];
45: VOLD = ksp->work[6];
46: WOLD = ksp->work[7];
47: WOOLD = ksp->work[8];
49: PCGetOperators(ksp->pc,&Amat,&Pmat);
51: ksp->its = 0;
53: VecSet(UOLD,0.0); /* u_old <- 0 */
54: VecSet(VOLD,0.0); /* v_old <- 0 */
55: VecSet(W,0.0); /* w <- 0 */
56: VecSet(WOLD,0.0); /* w_old <- 0 */
58: if (!ksp->guess_zero) {
59: KSP_MatMult(ksp,Amat,X,R); /* r <- b - A*x */
60: VecAYPX(R,-1.0,B);
61: } else {
62: VecCopy(B,R); /* r <- b (x is 0) */
63: }
64: KSP_PCApply(ksp,R,Z); /* z <- B*r */
65: VecNorm(Z,NORM_2,&np); /* np <- ||z|| */
66: VecDot(R,Z,&dp);
68: if (PetscRealPart(dp) < minres->haptol && np > minres->haptol) {
69: if (ksp->errorifnotconverged) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_CONV_FAILED,"Detected indefinite operator %g tolerance %g",(double)PetscRealPart(dp),(double)minres->haptol);
70: PetscInfo2(ksp,"Detected indefinite operator %g tolerance %g\n",(double)PetscRealPart(dp),(double)minres->haptol);
71: ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
72: return(0);
73: }
75: KSPLogResidualHistory(ksp,np);
76: KSPMonitor(ksp,0,np);
77: ksp->rnorm = np;
78: (*ksp->converged)(ksp,0,np,&ksp->reason,ksp->cnvP); /* test for convergence */
79: if (ksp->reason) return(0);
81: dp = PetscAbsScalar(dp);
82: dp = PetscSqrtScalar(dp);
83: beta = dp; /* beta <- sqrt(r'*z */
84: eta = beta;
86: VecCopy(R,V);
87: VecCopy(Z,U);
88: ibeta = 1.0 / beta;
89: VecScale(V,ibeta); /* v <- r / beta */
90: VecScale(U,ibeta); /* u <- z / beta */
92: i = 0;
93: do {
94: ksp->its = i+1;
96: /* Lanczos */
98: KSP_MatMult(ksp,Amat,U,R); /* r <- A*u */
99: VecDot(U,R,&alpha); /* alpha <- r'*u */
100: KSP_PCApply(ksp,R,Z); /* z <- B*r */
102: VecAXPY(R,-alpha,V); /* r <- r - alpha v */
103: VecAXPY(Z,-alpha,U); /* z <- z - alpha u */
104: VecAXPY(R,-beta,VOLD); /* r <- r - beta v_old */
105: VecAXPY(Z,-beta,UOLD); /* z <- z - beta u_old */
107: betaold = beta;
109: VecDot(R,Z,&dp);
110: dp = PetscAbsScalar(dp);
111: beta = PetscSqrtScalar(dp); /* beta <- sqrt(r'*z) */
113: /* QR factorisation */
115: coold = cold; cold = c; soold = sold; sold = s;
117: rho0 = cold * alpha - coold * sold * betaold;
118: rho1 = PetscSqrtScalar(rho0*rho0 + beta*beta);
119: rho2 = sold * alpha + coold * cold * betaold;
120: rho3 = soold * betaold;
122: /* Givens rotation */
124: c = rho0 / rho1;
125: s = beta / rho1;
127: /* Update */
129: VecCopy(WOLD,WOOLD); /* w_oold <- w_old */
130: VecCopy(W,WOLD); /* w_old <- w */
132: VecCopy(U,W); /* w <- u */
133: VecAXPY(W,-rho2,WOLD); /* w <- w - rho2 w_old */
134: VecAXPY(W,-rho3,WOOLD); /* w <- w - rho3 w_oold */
135: irho1 = 1.0 / rho1;
136: VecScale(W,irho1); /* w <- w / rho1 */
138: ceta = c * eta;
139: VecAXPY(X,ceta,W); /* x <- x + c eta w */
141: /*
142: when dp is really small we have either convergence or an indefinite operator so compute true
143: residual norm to check for convergence
144: */
145: if (PetscRealPart(dp) < minres->haptol) {
146: PetscInfo2(ksp,"Possible indefinite operator %g tolerance %g\n",(double)PetscRealPart(dp),(double)minres->haptol);
147: KSP_MatMult(ksp,Amat,X,VOLD);
148: VecAXPY(VOLD,none,B);
149: VecNorm(VOLD,NORM_2,&np);
150: } else {
151: /* otherwise compute new residual norm via recurrence relation */
152: np = ksp->rnorm * PetscAbsScalar(s);
153: }
155: ksp->rnorm = np;
156: KSPLogResidualHistory(ksp,np);
157: KSPMonitor(ksp,i+1,np);
158: (*ksp->converged)(ksp,i+1,np,&ksp->reason,ksp->cnvP); /* test for convergence */
159: if (ksp->reason) break;
161: if (PetscRealPart(dp) < minres->haptol) {
162: if (ksp->errorifnotconverged) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_CONV_FAILED,"Detected indefinite operator %g tolerance %g",(double)PetscRealPart(dp),(double)minres->haptol);
163: PetscInfo2(ksp,"Detected indefinite operator %g tolerance %g\n",(double)PetscRealPart(dp),(double)minres->haptol);
164: ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
165: break;
166: }
168: eta = -s * eta;
169: VecCopy(V,VOLD);
170: VecCopy(U,UOLD);
171: VecCopy(R,V);
172: VecCopy(Z,U);
173: ibeta = 1.0 / beta;
174: VecScale(V,ibeta); /* v <- r / beta */
175: VecScale(U,ibeta); /* u <- z / beta */
177: i++;
178: } while (i<ksp->max_it);
179: if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
180: return(0);
181: }
183: /*MC
184: KSPMINRES - This code implements the MINRES (Minimum Residual) method.
186: Options Database Keys:
187: . see KSPSolve()
189: Level: beginner
191: Notes:
192: The operator and the preconditioner must be symmetric and the preconditioner must
193: be positive definite for this method.
194: Supports only left preconditioning.
196: Reference: Paige & Saunders, 1975.
198: Contributed by: Robert Scheichl: maprs@maths.bath.ac.uk
200: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPCG, KSPCR
201: M*/
202: PETSC_EXTERN PetscErrorCode KSPCreate_MINRES(KSP ksp)
203: {
204: KSP_MINRES *minres;
208: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
209: PetscNewLog(ksp,&minres);
211: /* this parameter is arbitrary; but e-50 didn't work for __float128 in one example */
212: #if defined(PETSC_USE_REAL___FLOAT128)
213: minres->haptol = 1.e-100;
214: #elif defined(PETSC_USE_REAL_SINGLE)
215: minres->haptol = 1.e-25;
216: #else
217: minres->haptol = 1.e-50;
218: #endif
219: ksp->data = (void*)minres;
221: /*
222: Sets the functions that are associated with this data structure
223: (in C++ this is the same as defining virtual functions)
224: */
225: ksp->ops->setup = KSPSetUp_MINRES;
226: ksp->ops->solve = KSPSolve_MINRES;
227: ksp->ops->destroy = KSPDestroyDefault;
228: ksp->ops->setfromoptions = 0;
229: ksp->ops->buildsolution = KSPBuildSolutionDefault;
230: ksp->ops->buildresidual = KSPBuildResidualDefault;
231: return(0);
232: }