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