Actual source code: symmlq.c
petsc-3.5.4 2015-05-23
2: #include <petsc-private/kspimpl.h>
4: typedef struct {
5: PetscReal haptol;
6: } KSP_SYMMLQ;
10: PetscErrorCode KSPSetUp_SYMMLQ(KSP ksp)
11: {
15: KSPSetWorkVecs(ksp,9);
16: return(0);
17: }
21: PetscErrorCode KSPSolve_SYMMLQ(KSP ksp)
22: {
24: PetscInt i;
25: PetscScalar alpha,beta,ibeta,betaold,beta1,ceta = 0,ceta_oold = 0.0, ceta_old = 0.0,ceta_bar;
26: PetscScalar c = 1.0,cold=1.0,s=0.0,sold=0.0,coold,soold,rho0,rho1,rho2,rho3;
27: PetscScalar dp = 0.0;
28: PetscReal np,s_prod;
29: Vec X,B,R,Z,U,V,W,UOLD,VOLD,Wbar;
30: Mat Amat,Pmat;
31: KSP_SYMMLQ *symmlq = (KSP_SYMMLQ*)ksp->data;
32: PetscBool diagonalscale;
35: PCGetDiagonalScale(ksp->pc,&diagonalscale);
36: if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
38: X = ksp->vec_sol;
39: B = ksp->vec_rhs;
40: R = ksp->work[0];
41: Z = ksp->work[1];
42: U = ksp->work[2];
43: V = ksp->work[3];
44: W = ksp->work[4];
45: UOLD = ksp->work[5];
46: VOLD = ksp->work[6];
47: Wbar = ksp->work[7];
49: PCGetOperators(ksp->pc,&Amat,&Pmat);
51: ksp->its = 0;
53: VecSet(UOLD,0.0); /* u_old <- zeros; */
54: VecCopy(UOLD,VOLD); /* v_old <- u_old; */
55: VecCopy(UOLD,W); /* w <- u_old; */
56: VecCopy(UOLD,Wbar); /* w_bar <- u_old; */
57: if (!ksp->guess_zero) {
58: KSP_MatMult(ksp,Amat,X,R); /* r <- b - A*x */
59: VecAYPX(R,-1.0,B);
60: } else {
61: VecCopy(B,R); /* r <- b (x is 0) */
62: }
64: KSP_PCApply(ksp,R,Z); /* z <- B*r */
65: VecDot(R,Z,&dp); /* dp = r'*z; */
66: if (PetscAbsScalar(dp) < symmlq->haptol) {
67: PetscInfo2(ksp,"Detected happy breakdown %g tolerance %g\n",(double)PetscAbsScalar(dp),(double)symmlq->haptol);
68: ksp->rnorm = 0.0; /* what should we really put here? */
69: ksp->reason = KSP_CONVERGED_HAPPY_BREAKDOWN; /* bugfix proposed by Lourens (lourens.vanzanen@shell.com) */
70: return(0);
71: }
73: #if !defined(PETSC_USE_COMPLEX)
74: if (dp < 0.0) {
75: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
76: return(0);
77: }
78: #endif
79: dp = PetscSqrtScalar(dp);
80: beta = dp; /* beta <- sqrt(r'*z) */
81: beta1 = beta;
82: s_prod = PetscAbsScalar(beta1);
84: VecCopy(R,V); /* v <- r; */
85: VecCopy(Z,U); /* u <- z; */
86: ibeta = 1.0 / beta;
87: VecScale(V,ibeta); /* v <- ibeta*v; */
88: VecScale(U,ibeta); /* u <- ibeta*u; */
89: VecCopy(U,Wbar); /* w_bar <- u; */
90: VecNorm(Z,NORM_2,&np); /* np <- ||z|| */
91: KSPLogResidualHistory(ksp,np);
92: KSPMonitor(ksp,0,np);
93: ksp->rnorm = np;
94: (*ksp->converged)(ksp,0,np,&ksp->reason,ksp->cnvP); /* test for convergence */
95: if (ksp->reason) return(0);
97: i = 0; ceta = 0.;
98: do {
99: ksp->its = i+1;
101: /* Update */
102: if (ksp->its > 1) {
103: VecCopy(V,VOLD); /* v_old <- v; */
104: VecCopy(U,UOLD); /* u_old <- u; */
106: VecCopy(R,V);
107: VecScale(V,1.0/beta); /* v <- ibeta*r; */
108: VecCopy(Z,U);
109: VecScale(U,1.0/beta); /* u <- ibeta*z; */
111: VecCopy(Wbar,W);
112: VecScale(W,c);
113: VecAXPY(W,s,U); /* w <- c*w_bar + s*u; (w_k) */
114: VecScale(Wbar,-s);
115: VecAXPY(Wbar,c,U); /* w_bar <- -s*w_bar + c*u; (w_bar_(k+1)) */
116: VecAXPY(X,ceta,W); /* x <- x + ceta * w; (xL_k) */
118: ceta_oold = ceta_old;
119: ceta_old = ceta;
120: }
122: /* Lanczos */
123: KSP_MatMult(ksp,Amat,U,R); /* r <- Amat*u; */
124: VecDot(U,R,&alpha); /* alpha <- u'*r; */
125: KSP_PCApply(ksp,R,Z); /* z <- B*r; */
127: VecAXPY(R,-alpha,V); /* r <- r - alpha* v; */
128: VecAXPY(Z,-alpha,U); /* z <- z - alpha* u; */
129: VecAXPY(R,-beta,VOLD); /* r <- r - beta * v_old; */
130: VecAXPY(Z,-beta,UOLD); /* z <- z - beta * u_old; */
131: betaold = beta; /* beta_k */
132: VecDot(R,Z,&dp); /* dp <- r'*z; */
133: if (PetscAbsScalar(dp) < symmlq->haptol) {
134: PetscInfo2(ksp,"Detected happy breakdown %g tolerance %g\n",(double)PetscAbsScalar(dp),(double)symmlq->haptol);
135: dp = 0.0;
136: }
138: #if !defined(PETSC_USE_COMPLEX)
139: if (dp < 0.0) {
140: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
141: break;
142: }
143: #endif
144: beta = PetscSqrtScalar(dp); /* beta = sqrt(dp); */
146: /* QR factorization */
147: coold = cold; cold = c; soold = sold; sold = s;
148: rho0 = cold * alpha - coold * sold * betaold; /* gamma_bar */
149: rho1 = PetscSqrtScalar(rho0*rho0 + beta*beta); /* gamma */
150: rho2 = sold * alpha + coold * cold * betaold; /* delta */
151: rho3 = soold * betaold; /* epsilon */
153: /* Givens rotation: [c -s; s c] (different from the Reference!) */
154: c = rho0 / rho1; s = beta / rho1;
156: if (ksp->its==1) ceta = beta1/rho1;
157: else ceta = -(rho2*ceta_old + rho3*ceta_oold)/rho1;
159: s_prod = s_prod*PetscAbsScalar(s);
160: if (c == 0.0) np = s_prod*1.e16;
161: else np = s_prod/PetscAbsScalar(c); /* residual norm for xc_k (CGNORM) */
163: ksp->rnorm = np;
164: KSPLogResidualHistory(ksp,np);
165: KSPMonitor(ksp,i+1,np);
166: (*ksp->converged)(ksp,i+1,np,&ksp->reason,ksp->cnvP); /* test for convergence */
167: if (ksp->reason) break;
168: i++;
169: } while (i<ksp->max_it);
171: /* move to the CG point: xc_(k+1) */
172: if (c == 0.0) ceta_bar = ceta*1.e15;
173: else ceta_bar = ceta/c;
175: VecAXPY(X,ceta_bar,Wbar); /* x <- x + ceta_bar*w_bar */
177: if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
178: return(0);
179: }
181: /*MC
182: KSPSYMMLQ - This code implements the SYMMLQ method.
184: Options Database Keys:
185: . see KSPSolve()
187: Level: beginner
189: Notes: The operator and the preconditioner must be symmetric for this method. The
190: preconditioner must be POSITIVE-DEFINITE.
192: Supports only left preconditioning.
194: Reference: Paige & Saunders, 1975.
196: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP
197: M*/
200: PETSC_EXTERN PetscErrorCode KSPCreate_SYMMLQ(KSP ksp)
201: {
202: KSP_SYMMLQ *symmlq;
206: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
208: PetscNewLog(ksp,&symmlq);
209: symmlq->haptol = 1.e-18;
210: ksp->data = (void*)symmlq;
212: /*
213: Sets the functions that are associated with this data structure
214: (in C++ this is the same as defining virtual functions)
215: */
216: ksp->ops->setup = KSPSetUp_SYMMLQ;
217: ksp->ops->solve = KSPSolve_SYMMLQ;
218: ksp->ops->destroy = KSPDestroyDefault;
219: ksp->ops->setfromoptions = 0;
220: ksp->ops->buildsolution = KSPBuildSolutionDefault;
221: ksp->ops->buildresidual = KSPBuildResidualDefault;
222: return(0);
223: }