Actual source code: symmlq.c
petsc-3.14.3 2021-01-09
2: #include <petsc/private/kspimpl.h>
4: typedef struct {
5: PetscReal haptol;
6: } KSP_SYMMLQ;
8: PetscErrorCode KSPSetUp_SYMMLQ(KSP ksp)
9: {
13: KSPSetWorkVecs(ksp,9);
14: return(0);
15: }
17: PetscErrorCode KSPSolve_SYMMLQ(KSP ksp)
18: {
20: PetscInt i;
21: PetscScalar alpha,beta,ibeta,betaold,beta1,ceta = 0,ceta_oold = 0.0, ceta_old = 0.0,ceta_bar;
22: PetscScalar c = 1.0,cold=1.0,s=0.0,sold=0.0,coold,soold,rho0,rho1,rho2,rho3;
23: PetscScalar dp = 0.0;
24: PetscReal np = 0.0,s_prod;
25: Vec X,B,R,Z,U,V,W,UOLD,VOLD,Wbar;
26: Mat Amat,Pmat;
27: KSP_SYMMLQ *symmlq = (KSP_SYMMLQ*)ksp->data;
28: PetscBool diagonalscale;
31: PCGetDiagonalScale(ksp->pc,&diagonalscale);
32: if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
34: X = ksp->vec_sol;
35: B = ksp->vec_rhs;
36: R = ksp->work[0];
37: Z = ksp->work[1];
38: U = ksp->work[2];
39: V = ksp->work[3];
40: W = ksp->work[4];
41: UOLD = ksp->work[5];
42: VOLD = ksp->work[6];
43: Wbar = ksp->work[7];
45: PCGetOperators(ksp->pc,&Amat,&Pmat);
47: ksp->its = 0;
49: VecSet(UOLD,0.0); /* u_old <- zeros; */
50: VecCopy(UOLD,VOLD); /* v_old <- u_old; */
51: VecCopy(UOLD,W); /* w <- u_old; */
52: VecCopy(UOLD,Wbar); /* w_bar <- u_old; */
53: if (!ksp->guess_zero) {
54: KSP_MatMult(ksp,Amat,X,R); /* r <- b - A*x */
55: VecAYPX(R,-1.0,B);
56: } else {
57: VecCopy(B,R); /* r <- b (x is 0) */
58: }
60: KSP_PCApply(ksp,R,Z); /* z <- B*r */
61: VecDot(R,Z,&dp); /* dp = r'*z; */
62: KSPCheckDot(ksp,dp);
63: if (PetscAbsScalar(dp) < symmlq->haptol) {
64: PetscInfo2(ksp,"Detected happy breakdown %g tolerance %g\n",(double)PetscAbsScalar(dp),(double)symmlq->haptol);
65: ksp->rnorm = 0.0; /* what should we really put here? */
66: ksp->reason = KSP_CONVERGED_HAPPY_BREAKDOWN; /* bugfix proposed by Lourens (lourens.vanzanen@shell.com) */
67: return(0);
68: }
70: #if !defined(PETSC_USE_COMPLEX)
71: if (dp < 0.0) {
72: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
73: return(0);
74: }
75: #endif
76: dp = PetscSqrtScalar(dp);
77: beta = dp; /* beta <- sqrt(r'*z) */
78: beta1 = beta;
79: s_prod = PetscAbsScalar(beta1);
81: VecCopy(R,V); /* v <- r; */
82: VecCopy(Z,U); /* u <- z; */
83: ibeta = 1.0 / beta;
84: VecScale(V,ibeta); /* v <- ibeta*v; */
85: VecScale(U,ibeta); /* u <- ibeta*u; */
86: VecCopy(U,Wbar); /* w_bar <- u; */
87: if (ksp->normtype != KSP_NORM_NONE) {
88: VecNorm(Z,NORM_2,&np); /* np <- ||z|| */
89: KSPCheckNorm(ksp,np);
90: }
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: KSPCheckDot(ksp,dp);
134: if (PetscAbsScalar(dp) < symmlq->haptol) {
135: PetscInfo2(ksp,"Detected happy breakdown %g tolerance %g\n",(double)PetscAbsScalar(dp),(double)symmlq->haptol);
136: dp = 0.0;
137: }
139: #if !defined(PETSC_USE_COMPLEX)
140: if (dp < 0.0) {
141: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
142: break;
143: }
144: #endif
145: beta = PetscSqrtScalar(dp); /* beta = sqrt(dp); */
147: /* QR factorization */
148: coold = cold; cold = c; soold = sold; sold = s;
149: rho0 = cold * alpha - coold * sold * betaold; /* gamma_bar */
150: rho1 = PetscSqrtScalar(rho0*rho0 + beta*beta); /* gamma */
151: rho2 = sold * alpha + coold * cold * betaold; /* delta */
152: rho3 = soold * betaold; /* epsilon */
154: /* Givens rotation: [c -s; s c] (different from the Reference!) */
155: c = rho0 / rho1; s = beta / rho1;
157: if (ksp->its==1) ceta = beta1/rho1;
158: else ceta = -(rho2*ceta_old + rho3*ceta_oold)/rho1;
160: s_prod = s_prod*PetscAbsScalar(s);
161: if (c == 0.0) np = s_prod*1.e16;
162: else np = s_prod/PetscAbsScalar(c); /* residual norm for xc_k (CGNORM) */
164: if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = np;
165: else ksp->rnorm = 0.0;
166: KSPLogResidualHistory(ksp,ksp->rnorm);
167: KSPMonitor(ksp,i+1,ksp->rnorm);
168: (*ksp->converged)(ksp,i+1,ksp->rnorm,&ksp->reason,ksp->cnvP); /* test for convergence */
169: if (ksp->reason) break;
170: i++;
171: } while (i<ksp->max_it);
173: /* move to the CG point: xc_(k+1) */
174: if (c == 0.0) ceta_bar = ceta*1.e15;
175: else ceta_bar = ceta/c;
177: VecAXPY(X,ceta_bar,Wbar); /* x <- x + ceta_bar*w_bar */
179: if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
180: return(0);
181: }
183: /*MC
184: KSPSYMMLQ - This code implements the SYMMLQ method.
186: Options Database Keys:
187: . see KSPSolve()
189: Level: beginner
191: Notes:
192: The operator and the preconditioner must be symmetric for this method. The
193: preconditioner must be POSITIVE-DEFINITE.
195: Supports only left preconditioning.
197: Reference: Paige & Saunders, 1975.
199: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP
200: M*/
201: PETSC_EXTERN PetscErrorCode KSPCreate_SYMMLQ(KSP ksp)
202: {
203: KSP_SYMMLQ *symmlq;
207: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
208: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);
210: PetscNewLog(ksp,&symmlq);
211: symmlq->haptol = 1.e-18;
212: ksp->data = (void*)symmlq;
214: /*
215: Sets the functions that are associated with this data structure
216: (in C++ this is the same as defining virtual functions)
217: */
218: ksp->ops->setup = KSPSetUp_SYMMLQ;
219: ksp->ops->solve = KSPSolve_SYMMLQ;
220: ksp->ops->destroy = KSPDestroyDefault;
221: ksp->ops->setfromoptions = NULL;
222: ksp->ops->buildsolution = KSPBuildSolutionDefault;
223: ksp->ops->buildresidual = KSPBuildResidualDefault;
224: return(0);
225: }