Actual source code: bcgs.c
petsc-3.13.6 2020-09-29
2: #include <../src/ksp/ksp/impls/bcgs/bcgsimpl.h>
4: PetscErrorCode KSPSetFromOptions_BCGS(PetscOptionItems *PetscOptionsObject,KSP ksp)
5: {
9: PetscOptionsHead(PetscOptionsObject,"KSP BCGS Options");
10: PetscOptionsTail();
11: return(0);
12: }
14: PetscErrorCode KSPSetUp_BCGS(KSP ksp)
15: {
19: KSPSetWorkVecs(ksp,6);
20: return(0);
21: }
24: PetscErrorCode KSPSolve_BCGS(KSP ksp)
25: {
27: PetscInt i;
28: PetscScalar rho,rhoold,alpha,beta,omega,omegaold,d1;
29: Vec X,B,V,P,R,RP,T,S;
30: PetscReal dp = 0.0,d2;
31: KSP_BCGS *bcgs = (KSP_BCGS*)ksp->data;
34: X = ksp->vec_sol;
35: B = ksp->vec_rhs;
36: R = ksp->work[0];
37: RP = ksp->work[1];
38: V = ksp->work[2];
39: T = ksp->work[3];
40: S = ksp->work[4];
41: P = ksp->work[5];
43: /* Compute initial preconditioned residual */
44: KSPInitialResidual(ksp,X,V,T,R,B);
46: /* with right preconditioning need to save initial guess to add to final solution */
47: if (ksp->pc_side == PC_RIGHT && !ksp->guess_zero) {
48: if (!bcgs->guess) {
49: VecDuplicate(X,&bcgs->guess);
50: }
51: VecCopy(X,bcgs->guess);
52: VecSet(X,0.0);
53: }
55: /* Test for nothing to do */
56: if (ksp->normtype != KSP_NORM_NONE) {
57: VecNorm(R,NORM_2,&dp);
58: KSPCheckNorm(ksp,dp);
59: }
60: PetscObjectSAWsTakeAccess((PetscObject)ksp);
61: ksp->its = 0;
62: ksp->rnorm = dp;
63: PetscObjectSAWsGrantAccess((PetscObject)ksp);
64: KSPLogResidualHistory(ksp,dp);
65: KSPMonitor(ksp,0,dp);
66: (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
67: if (ksp->reason) {
68: if (bcgs->guess) {
69: VecAXPY(X,1.0,bcgs->guess);
70: }
71: return(0);
72: }
74: /* Make the initial Rp == R */
75: VecCopy(R,RP);
77: rhoold = 1.0;
78: alpha = 1.0;
79: omegaold = 1.0;
80: VecSet(P,0.0);
81: VecSet(V,0.0);
83: i=0;
84: do {
85: VecDot(R,RP,&rho); /* rho <- (r,rp) */
86: beta = (rho/rhoold) * (alpha/omegaold);
87: VecAXPBYPCZ(P,1.0,-omegaold*beta,beta,R,V); /* p <- r - omega * beta* v + beta * p */
88: KSP_PCApplyBAorAB(ksp,P,V,T); /* v <- K p */
89: VecDot(V,RP,&d1);
90: KSPCheckDot(ksp,d1);
91: if (d1 == 0.0) {
92: if (ksp->errorifnotconverged) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged due to Nan or Inf inner product");
93: else {
94: ksp->reason = KSP_DIVERGED_NANORINF;
95: break;
96: }
97: }
98: alpha = rho / d1; /* a <- rho / (v,rp) */
99: VecWAXPY(S,-alpha,V,R); /* s <- r - a v */
100: KSP_PCApplyBAorAB(ksp,S,T,R); /* t <- K s */
101: VecDotNorm2(S,T,&d1,&d2);
102: if (d2 == 0.0) {
103: /* t is 0. if s is 0, then alpha v == r, and hence alpha p
104: may be our solution. Give it a try? */
105: VecDot(S,S,&d1);
106: if (d1 != 0.0) {
107: ksp->reason = KSP_DIVERGED_BREAKDOWN;
108: break;
109: }
110: VecAXPY(X,alpha,P); /* x <- x + a p */
111: PetscObjectSAWsTakeAccess((PetscObject)ksp);
112: ksp->its++;
113: ksp->rnorm = 0.0;
114: ksp->reason = KSP_CONVERGED_RTOL;
115: PetscObjectSAWsGrantAccess((PetscObject)ksp);
116: KSPLogResidualHistory(ksp,dp);
117: KSPMonitor(ksp,i+1,0.0);
118: break;
119: }
120: omega = d1 / d2; /* w <- (t's) / (t't) */
121: VecAXPBYPCZ(X,alpha,omega,1.0,P,S); /* x <- alpha * p + omega * s + x */
122: VecWAXPY(R,-omega,T,S); /* r <- s - w t */
123: if (ksp->normtype != KSP_NORM_NONE && ksp->chknorm < i+2) {
124: VecNorm(R,NORM_2,&dp);
125: KSPCheckNorm(ksp,dp);
126: }
128: rhoold = rho;
129: omegaold = omega;
131: PetscObjectSAWsTakeAccess((PetscObject)ksp);
132: ksp->its++;
133: ksp->rnorm = dp;
134: PetscObjectSAWsGrantAccess((PetscObject)ksp);
135: KSPLogResidualHistory(ksp,dp);
136: KSPMonitor(ksp,i+1,dp);
137: (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
138: if (ksp->reason) break;
139: if (rho == 0.0) {
140: ksp->reason = KSP_DIVERGED_BREAKDOWN;
141: break;
142: }
143: i++;
144: } while (i<ksp->max_it);
146: if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
148: KSPUnwindPreconditioner(ksp,X,T);
149: if (bcgs->guess) {
150: VecAXPY(X,1.0,bcgs->guess);
151: }
152: return(0);
153: }
155: PetscErrorCode KSPBuildSolution_BCGS(KSP ksp,Vec v,Vec *V)
156: {
158: KSP_BCGS *bcgs = (KSP_BCGS*)ksp->data;
161: if (ksp->pc_side == PC_RIGHT) {
162: if (v) {
163: KSP_PCApply(ksp,ksp->vec_sol,v);
164: if (bcgs->guess) {
165: VecAXPY(v,1.0,bcgs->guess);
166: }
167: *V = v;
168: } else SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Not working with right preconditioner");
169: } else {
170: if (v) {
171: VecCopy(ksp->vec_sol,v); *V = v;
172: } else *V = ksp->vec_sol;
173: }
174: return(0);
175: }
177: PetscErrorCode KSPReset_BCGS(KSP ksp)
178: {
179: KSP_BCGS *cg = (KSP_BCGS*)ksp->data;
183: VecDestroy(&cg->guess);
184: return(0);
185: }
187: PetscErrorCode KSPDestroy_BCGS(KSP ksp)
188: {
192: KSPReset_BCGS(ksp);
193: KSPDestroyDefault(ksp);
194: return(0);
195: }
197: /*MC
198: KSPBCGS - Implements the BiCGStab (Stabilized version of BiConjugate Gradient) method.
200: Options Database Keys:
201: . see KSPSolve()
203: Level: beginner
205: Notes:
206: See KSPBCGSL for additional stabilization
207: Supports left and right preconditioning but not symmetric
209: References:
210: . 1. - van der Vorst, SIAM J. Sci. Stat. Comput., 1992.
212: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPBICG, KSPBCGSL, KSPFBICG, KSPSetPCSide()
213: M*/
214: PETSC_EXTERN PetscErrorCode KSPCreate_BCGS(KSP ksp)
215: {
217: KSP_BCGS *bcgs;
220: PetscNewLog(ksp,&bcgs);
222: ksp->data = bcgs;
223: ksp->ops->setup = KSPSetUp_BCGS;
224: ksp->ops->solve = KSPSolve_BCGS;
225: ksp->ops->destroy = KSPDestroy_BCGS;
226: ksp->ops->reset = KSPReset_BCGS;
227: ksp->ops->buildsolution = KSPBuildSolution_BCGS;
228: ksp->ops->buildresidual = KSPBuildResidualDefault;
229: ksp->ops->setfromoptions = KSPSetFromOptions_BCGS;
231: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
232: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);
233: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);
234: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,1);
235: return(0);
236: }