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: dp = 0.0;
76: }
78: #if !defined(PETSC_USE_COMPLEX)
79: if (dp < 0.0) {
80: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
81: return(0);
82: }
83: #endif
84: dp = PetscSqrtScalar(dp);
85: beta = dp; /* beta <- sqrt(r'*z) */
86: beta1 = beta;
87: s_prod = PetscAbsScalar(beta1);
89: VecCopy(R,V); /* v <- r; */
90: VecCopy(Z,U); /* u <- z; */
91: ibeta = 1.0 / beta;
92: VecScale(V,ibeta); /* v <- ibeta*v; */
93: VecScale(U,ibeta); /* u <- ibeta*u; */
94: VecCopy(U,Wbar); /* w_bar <- u; */
95: VecNorm(Z,NORM_2,&np); /* np <- ||z|| */
96: KSPLogResidualHistory(ksp,np);
97: KSPMonitor(ksp,0,np); /* call any registered monitor routines */
98: ksp->rnorm = np;
99: (*ksp->converged)(ksp,0,np,&ksp->reason,ksp->cnvP); /* test for convergence */
100: if (ksp->reason) return(0);
102: i = 0; ceta = 0.;
103: do {
104: ksp->its = i+1;
106: /* Update */
107: if (ksp->its > 1){
108: VecCopy(V,VOLD); /* v_old <- v; */
109: VecCopy(U,UOLD); /* u_old <- u; */
110:
111: VecCopy(R,V);
112: VecScale(V,1.0/beta); /* v <- ibeta*r; */
113: VecCopy(Z,U);
114: VecScale(U,1.0/beta); /* u <- ibeta*z; */
116: VecCopy(Wbar,W);
117: VecScale(W,c);
118: VecAXPY(W,s,U); /* w <- c*w_bar + s*u; (w_k) */
119: VecScale(Wbar,-s);
120: VecAXPY(Wbar,c,U); /* w_bar <- -s*w_bar + c*u; (w_bar_(k+1)) */
121: VecAXPY(X,ceta,W); /* x <- x + ceta * w; (xL_k) */
123: ceta_oold = ceta_old;
124: ceta_old = ceta;
125: }
127: /* Lanczos */
128: KSP_MatMult(ksp,Amat,U,R); /* r <- Amat*u; */
129: VecDot(U,R,&alpha); /* alpha <- u'*r; */
130: KSP_PCApply(ksp,R,Z); /* z <- B*r; */
132: VecAXPY(R,-alpha,V); /* r <- r - alpha* v; */
133: VecAXPY(Z,-alpha,U); /* z <- z - alpha* u; */
134: VecAXPY(R,-beta,VOLD); /* r <- r - beta * v_old; */
135: VecAXPY(Z,-beta,UOLD); /* z <- z - beta * u_old; */
136: betaold = beta; /* beta_k */
137: VecDot(R,Z,&dp); /* dp <- r'*z; */
138: if (PetscAbsScalar(dp) < symmlq->haptol) {
139: PetscInfo2(ksp,"Detected happy breakdown %G tolerance %G\n",PetscAbsScalar(dp),symmlq->haptol);
140: dp = 0.0;
141: }
143: #if !defined(PETSC_USE_COMPLEX)
144: if (dp < 0.0) {
145: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
146: break;
147: }
148: #endif
149: beta = PetscSqrtScalar(dp); /* beta = sqrt(dp); */
151: /* QR factorization */
152: coold = cold; cold = c; soold = sold; sold = s;
153: rho0 = cold * alpha - coold * sold * betaold; /* gamma_bar */
154: rho1 = PetscSqrtScalar(rho0*rho0 + beta*beta); /* gamma */
155: rho2 = sold * alpha + coold * cold * betaold; /* delta */
156: rho3 = soold * betaold; /* epsilon */
158: /* Givens rotation: [c -s; s c] (different from the Reference!) */
159: c = rho0 / rho1; s = beta / rho1;
161: if (ksp->its==1){
162: ceta = beta1/rho1;
163: } else {
164: ceta = -(rho2*ceta_old + rho3*ceta_oold)/rho1;
165: }
166:
167: s_prod = s_prod*PetscAbsScalar(s);
168: if (c == 0.0){
169: np = s_prod*1.e16;
170: } else {
171: np = s_prod/PetscAbsScalar(c); /* residual norm for xc_k (CGNORM) */
172: }
173: ksp->rnorm = np;
174: KSPLogResidualHistory(ksp,np);
175: KSPMonitor(ksp,i+1,np);
176: (*ksp->converged)(ksp,i+1,np,&ksp->reason,ksp->cnvP); /* test for convergence */
177: if (ksp->reason) break;
178: i++;
179: } while (i<ksp->max_it);
181: /* move to the CG point: xc_(k+1) */
182: if (c == 0.0){
183: ceta_bar = ceta*1.e15;
184: } else {
185: ceta_bar = ceta/c;
186: }
187: VecAXPY(X,ceta_bar,Wbar); /* x <- x + ceta_bar*w_bar */
189: if (i >= ksp->max_it) {
190: ksp->reason = KSP_DIVERGED_ITS;
191: }
192: return(0);
193: }
195: /*MC
196: KSPSYMMLQ - This code implements the SYMMLQ method.
198: Options Database Keys:
199: . see KSPSolve()
201: Level: beginner
203: Notes: Reference: Paige & Saunders, 1975.
205: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP
206: M*/
210: PetscErrorCode PETSCKSP_DLLEXPORT KSPCreate_SYMMLQ(KSP ksp)
211: {
212: KSP_SYMMLQ *symmlq;
216: ksp->pc_side = PC_LEFT;
218: PetscNew(KSP_SYMMLQ,&symmlq);
219: symmlq->haptol = 1.e-18;
220: ksp->data = (void*)symmlq;
222: /*
223: Sets the functions that are associated with this data structure
224: (in C++ this is the same as defining virtual functions)
225: */
226: ksp->ops->setup = KSPSetUp_SYMMLQ;
227: ksp->ops->solve = KSPSolve_SYMMLQ;
228: ksp->ops->destroy = KSPDefaultDestroy;
229: ksp->ops->setfromoptions = 0;
230: ksp->ops->buildsolution = KSPDefaultBuildSolution;
231: ksp->ops->buildresidual = KSPDefaultBuildResidual;
232: return(0);
233: }