Actual source code: bcgs.c
petsc-3.10.5 2019-03-28
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: }
59: PetscObjectSAWsTakeAccess((PetscObject)ksp);
60: ksp->its = 0;
61: ksp->rnorm = dp;
62: PetscObjectSAWsGrantAccess((PetscObject)ksp);
63: KSPLogResidualHistory(ksp,dp);
64: KSPMonitor(ksp,0,dp);
65: (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
66: if (ksp->reason) {
67: if (bcgs->guess) {
68: VecAXPY(X,1.0,bcgs->guess);
69: }
70: return(0);
71: }
73: /* Make the initial Rp == R */
74: VecCopy(R,RP);
76: rhoold = 1.0;
77: alpha = 1.0;
78: omegaold = 1.0;
79: VecSet(P,0.0);
80: VecSet(V,0.0);
82: i=0;
83: do {
84: VecDot(R,RP,&rho); /* rho <- (r,rp) */
85: beta = (rho/rhoold) * (alpha/omegaold);
86: VecAXPBYPCZ(P,1.0,-omegaold*beta,beta,R,V); /* p <- r - omega * beta* v + beta * p */
87: KSP_PCApplyBAorAB(ksp,P,V,T); /* v <- K p */
88: VecDot(V,RP,&d1);
89: if (d1 == 0.0) {
90: if (ksp->errorifnotconverged) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged due to Nan or Inf inner product");
91: else {
92: ksp->reason = KSP_DIVERGED_NANORINF;
93: break;
94: }
95: }
96: alpha = rho / d1; /* a <- rho / (v,rp) */
97: VecWAXPY(S,-alpha,V,R); /* s <- r - a v */
98: KSP_PCApplyBAorAB(ksp,S,T,R); /* t <- K s */
99: VecDotNorm2(S,T,&d1,&d2);
100: if (d2 == 0.0) {
101: /* t is 0. if s is 0, then alpha v == r, and hence alpha p
102: may be our solution. Give it a try? */
103: VecDot(S,S,&d1);
104: if (d1 != 0.0) {
105: ksp->reason = KSP_DIVERGED_BREAKDOWN;
106: break;
107: }
108: VecAXPY(X,alpha,P); /* x <- x + a p */
109: PetscObjectSAWsTakeAccess((PetscObject)ksp);
110: ksp->its++;
111: ksp->rnorm = 0.0;
112: ksp->reason = KSP_CONVERGED_RTOL;
113: PetscObjectSAWsGrantAccess((PetscObject)ksp);
114: KSPLogResidualHistory(ksp,dp);
115: KSPMonitor(ksp,i+1,0.0);
116: break;
117: }
118: omega = d1 / d2; /* w <- (t's) / (t't) */
119: VecAXPBYPCZ(X,alpha,omega,1.0,P,S); /* x <- alpha * p + omega * s + x */
120: VecWAXPY(R,-omega,T,S); /* r <- s - w t */
121: if (ksp->normtype != KSP_NORM_NONE && ksp->chknorm < i+2) {
122: VecNorm(R,NORM_2,&dp);
123: }
125: rhoold = rho;
126: omegaold = omega;
128: PetscObjectSAWsTakeAccess((PetscObject)ksp);
129: ksp->its++;
130: ksp->rnorm = dp;
131: PetscObjectSAWsGrantAccess((PetscObject)ksp);
132: KSPLogResidualHistory(ksp,dp);
133: KSPMonitor(ksp,i+1,dp);
134: (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
135: if (ksp->reason) break;
136: if (rho == 0.0) {
137: ksp->reason = KSP_DIVERGED_BREAKDOWN;
138: break;
139: }
140: i++;
141: } while (i<ksp->max_it);
143: if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
145: KSPUnwindPreconditioner(ksp,X,T);
146: if (bcgs->guess) {
147: VecAXPY(X,1.0,bcgs->guess);
148: }
149: return(0);
150: }
152: PetscErrorCode KSPBuildSolution_BCGS(KSP ksp,Vec v,Vec *V)
153: {
155: KSP_BCGS *bcgs = (KSP_BCGS*)ksp->data;
158: if (ksp->pc_side == PC_RIGHT) {
159: if (v) {
160: KSP_PCApply(ksp,ksp->vec_sol,v);
161: if (bcgs->guess) {
162: VecAXPY(v,1.0,bcgs->guess);
163: }
164: *V = v;
165: } else SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Not working with right preconditioner");
166: } else {
167: if (v) {
168: VecCopy(ksp->vec_sol,v); *V = v;
169: } else *V = ksp->vec_sol;
170: }
171: return(0);
172: }
174: PetscErrorCode KSPReset_BCGS(KSP ksp)
175: {
176: KSP_BCGS *cg = (KSP_BCGS*)ksp->data;
180: VecDestroy(&cg->guess);
181: return(0);
182: }
184: PetscErrorCode KSPDestroy_BCGS(KSP ksp)
185: {
189: KSPReset_BCGS(ksp);
190: KSPDestroyDefault(ksp);
191: return(0);
192: }
194: /*MC
195: KSPBCGS - Implements the BiCGStab (Stabilized version of BiConjugate Gradient) method.
197: Options Database Keys:
198: . see KSPSolve()
200: Level: beginner
202: Notes:
203: See KSPBCGSL for additional stabilization
204: Supports left and right preconditioning but not symmetric
206: References:
207: . 1. - van der Vorst, SIAM J. Sci. Stat. Comput., 1992.
209: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPBICG, KSPBCGSL, KSPFBICG, KSPSetPCSide()
210: M*/
211: PETSC_EXTERN PetscErrorCode KSPCreate_BCGS(KSP ksp)
212: {
214: KSP_BCGS *bcgs;
217: PetscNewLog(ksp,&bcgs);
219: ksp->data = bcgs;
220: ksp->ops->setup = KSPSetUp_BCGS;
221: ksp->ops->solve = KSPSolve_BCGS;
222: ksp->ops->destroy = KSPDestroy_BCGS;
223: ksp->ops->reset = KSPReset_BCGS;
224: ksp->ops->buildsolution = KSPBuildSolution_BCGS;
225: ksp->ops->buildresidual = KSPBuildResidualDefault;
226: ksp->ops->setfromoptions = KSPSetFromOptions_BCGS;
228: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
229: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);
230: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);
231: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,1);
232: return(0);
233: }