Actual source code: symmlq.c
1: #define PETSCKSP_DLL
3: #include src/ksp/ksp/kspimpl.h
5: typedef struct {
6: PetscReal haptol;
7: } KSP_SYMMLQ;
11: PetscErrorCode KSPSetUp_SYMMLQ(KSP ksp)
12: {
16: if (ksp->pc_side == PC_RIGHT) {
17: SETERRQ(PETSC_ERR_SUP,"No right preconditioning for KSPSYMMLQ");
18: } else if (ksp->pc_side == PC_SYMMETRIC) {
19: SETERRQ(PETSC_ERR_SUP,"No symmetric preconditioning for KSPSYMMLQ");
20: }
21: KSPDefaultGetWork(ksp,9);
22: return(0);
23: }
27: PetscErrorCode KSPSolve_SYMMLQ(KSP ksp)
28: {
30: PetscInt i;
31: PetscScalar alpha,beta,ibeta,betaold,beta1,ceta = 0,ceta_oold = 0.0, ceta_old = 0.0,ceta_bar;
32: PetscScalar c=1.0,cold=1.0,s=0.0,sold=0.0,coold,soold,rho0,rho1,rho2,rho3;
33: PetscScalar dp = 0.0;
34: PetscReal np,s_prod;
35: Vec X,B,R,Z,U,V,W,UOLD,VOLD,Wbar;
36: Mat Amat,Pmat;
37: MatStructure pflag;
38: KSP_SYMMLQ *symmlq = (KSP_SYMMLQ*)ksp->data;
39: PetscTruth diagonalscale;
42: PCDiagonalScale(ksp->pc,&diagonalscale);
43: if (diagonalscale) SETERRQ1(PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",ksp->type_name);
45: X = ksp->vec_sol;
46: B = ksp->vec_rhs;
47: R = ksp->work[0];
48: Z = ksp->work[1];
49: U = ksp->work[2];
50: V = ksp->work[3];
51: W = ksp->work[4];
52: UOLD = ksp->work[5];
53: VOLD = ksp->work[6];
54: Wbar = ksp->work[7];
55:
56: PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);
58: ksp->its = 0;
60: VecSet(UOLD,0.0); /* u_old <- zeros; */
61: VecCopy(UOLD,VOLD); /* v_old <- u_old; */
62: VecCopy(UOLD,W); /* w <- u_old; */
63: VecCopy(UOLD,Wbar); /* w_bar <- u_old; */
64: if (!ksp->guess_zero) {
65: KSP_MatMult(ksp,Amat,X,R); /* r <- b - A*x */
66: VecAYPX(R,-1.0,B);
67: } else {
68: VecCopy(B,R); /* r <- b (x is 0) */
69: }
71: KSP_PCApply(ksp,R,Z); /* z <- B*r */
72: VecDot(R,Z,&dp); /* dp = r'*z; */
73: if (PetscAbsScalar(dp) < symmlq->haptol) {
74: PetscInfo2(ksp,"Detected happy breakdown %G tolerance %G\n",PetscAbsScalar(dp),symmlq->haptol);
75: ksp->rnorm = 0.0; /* what should we really put here? */
76: ksp->reason = KSP_CONVERGED_HAPPY_BREAKDOWN; /* bugfix proposed by Lourens (lourens.vanzanen@shell.com) */
77: return(0);
78: }
80: #if !defined(PETSC_USE_COMPLEX)
81: if (dp < 0.0) {
82: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
83: return(0);
84: }
85: #endif
86: dp = PetscSqrtScalar(dp);
87: beta = dp; /* beta <- sqrt(r'*z) */
88: beta1 = beta;
89: s_prod = PetscAbsScalar(beta1);
91: VecCopy(R,V); /* v <- r; */
92: VecCopy(Z,U); /* u <- z; */
93: ibeta = 1.0 / beta;
94: VecScale(V,ibeta); /* v <- ibeta*v; */
95: VecScale(U,ibeta); /* u <- ibeta*u; */
96: VecCopy(U,Wbar); /* w_bar <- u; */
97: VecNorm(Z,NORM_2,&np); /* np <- ||z|| */
98: KSPLogResidualHistory(ksp,np);
99: KSPMonitor(ksp,0,np); /* call any registered monitor routines */
100: ksp->rnorm = np;
101: (*ksp->converged)(ksp,0,np,&ksp->reason,ksp->cnvP); /* test for convergence */
102: if (ksp->reason) return(0);
104: i = 0; ceta = 0.;
105: do {
106: ksp->its = i+1;
108: /* Update */
109: if (ksp->its > 1){
110: VecCopy(V,VOLD); /* v_old <- v; */
111: VecCopy(U,UOLD); /* u_old <- u; */
112:
113: VecCopy(R,V);
114: VecScale(V,1.0/beta); /* v <- ibeta*r; */
115: VecCopy(Z,U);
116: VecScale(U,1.0/beta); /* u <- ibeta*z; */
118: VecCopy(Wbar,W);
119: VecScale(W,c);
120: VecAXPY(W,s,U); /* w <- c*w_bar + s*u; (w_k) */
121: VecScale(Wbar,-s);
122: VecAXPY(Wbar,c,U); /* w_bar <- -s*w_bar + c*u; (w_bar_(k+1)) */
123: VecAXPY(X,ceta,W); /* x <- x + ceta * w; (xL_k) */
125: ceta_oold = ceta_old;
126: ceta_old = ceta;
127: }
129: /* Lanczos */
130: KSP_MatMult(ksp,Amat,U,R); /* r <- Amat*u; */
131: VecDot(U,R,&alpha); /* alpha <- u'*r; */
132: KSP_PCApply(ksp,R,Z); /* z <- B*r; */
134: VecAXPY(R,-alpha,V); /* r <- r - alpha* v; */
135: VecAXPY(Z,-alpha,U); /* z <- z - alpha* u; */
136: VecAXPY(R,-beta,VOLD); /* r <- r - beta * v_old; */
137: VecAXPY(Z,-beta,UOLD); /* z <- z - beta * u_old; */
138: betaold = beta; /* beta_k */
139: VecDot(R,Z,&dp); /* dp <- r'*z; */
140: if (PetscAbsScalar(dp) < symmlq->haptol) {
141: PetscInfo2(ksp,"Detected happy breakdown %G tolerance %G\n",PetscAbsScalar(dp),symmlq->haptol);
142: dp = 0.0;
143: }
145: #if !defined(PETSC_USE_COMPLEX)
146: if (dp < 0.0) {
147: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
148: break;
149: }
150: #endif
151: beta = PetscSqrtScalar(dp); /* beta = sqrt(dp); */
153: /* QR factorization */
154: coold = cold; cold = c; soold = sold; sold = s;
155: rho0 = cold * alpha - coold * sold * betaold; /* gamma_bar */
156: rho1 = PetscSqrtScalar(rho0*rho0 + beta*beta); /* gamma */
157: rho2 = sold * alpha + coold * cold * betaold; /* delta */
158: rho3 = soold * betaold; /* epsilon */
160: /* Givens rotation: [c -s; s c] (different from the Reference!) */
161: c = rho0 / rho1; s = beta / rho1;
163: if (ksp->its==1){
164: ceta = beta1/rho1;
165: } else {
166: ceta = -(rho2*ceta_old + rho3*ceta_oold)/rho1;
167: }
168:
169: s_prod = s_prod*PetscAbsScalar(s);
170: if (c == 0.0){
171: np = s_prod*1.e16;
172: } else {
173: np = s_prod/PetscAbsScalar(c); /* residual norm for xc_k (CGNORM) */
174: }
175: ksp->rnorm = np;
176: KSPLogResidualHistory(ksp,np);
177: KSPMonitor(ksp,i+1,np);
178: (*ksp->converged)(ksp,i+1,np,&ksp->reason,ksp->cnvP); /* test for convergence */
179: if (ksp->reason) break;
180: i++;
181: } while (i<ksp->max_it);
183: /* move to the CG point: xc_(k+1) */
184: if (c == 0.0){
185: ceta_bar = ceta*1.e15;
186: } else {
187: ceta_bar = ceta/c;
188: }
189: VecAXPY(X,ceta_bar,Wbar); /* x <- x + ceta_bar*w_bar */
191: if (i >= ksp->max_it) {
192: ksp->reason = KSP_DIVERGED_ITS;
193: }
194: return(0);
195: }
197: /*MC
198: KSPSYMMLQ - This code implements the SYMMLQ method.
200: Options Database Keys:
201: . see KSPSolve()
203: Level: beginner
205: Notes: Reference: Paige & Saunders, 1975.
207: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP
208: M*/
212: PetscErrorCode KSPCreate_SYMMLQ(KSP ksp)
213: {
214: KSP_SYMMLQ *symmlq;
218: ksp->pc_side = PC_LEFT;
220: PetscNew(KSP_SYMMLQ,&symmlq);
221: symmlq->haptol = 1.e-18;
222: ksp->data = (void*)symmlq;
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_SYMMLQ;
229: ksp->ops->solve = KSPSolve_SYMMLQ;
230: ksp->ops->destroy = KSPDefaultDestroy;
231: ksp->ops->setfromoptions = 0;
232: ksp->ops->buildsolution = KSPDefaultBuildSolution;
233: ksp->ops->buildresidual = KSPDefaultBuildResidual;
234: return(0);
235: }