Actual source code: bcgs.c
1: #define PETSCKSP_DLL
3: #include private/kspimpl.h
7: static PetscErrorCode KSPSetUp_BCGS(KSP ksp)
8: {
12: if (ksp->pc_side == PC_SYMMETRIC) {
13: SETERRQ(PETSC_ERR_SUP,"no symmetric preconditioning for KSPBCGS");
14: }
15: KSPDefaultGetWork(ksp,6);
16: return(0);
17: }
21: static PetscErrorCode KSPSolve_BCGS(KSP ksp)
22: {
24: PetscInt i;
25: PetscScalar rho,rhoold,alpha,beta,omega,omegaold,d1,d2;
26: Vec X,B,V,P,R,RP,T,S;
27: PetscReal dp = 0.0;
30: if (ksp->normtype == KSP_NORM_NATURAL) SETERRQ(PETSC_ERR_SUP,"Cannot use natural residual norm with KSPBCGS");
31: if (ksp->normtype == KSP_NORM_PRECONDITIONED && ksp->pc_side != PC_LEFT) SETERRQ(PETSC_ERR_SUP,"Use -ksp_norm_type unpreconditioned for right preconditioning and KSPBCGS");
32: if (ksp->normtype == KSP_NORM_UNPRECONDITIONED && ksp->pc_side != PC_RIGHT) SETERRQ(PETSC_ERR_SUP,"Use -ksp_norm_type preconditioned for left preconditioning and KSPBCGS");
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: /* Test for nothing to do */
47: if (ksp->normtype != KSP_NORM_NO) {
48: VecNorm(R,NORM_2,&dp);
49: }
50: PetscObjectTakeAccess(ksp);
51: ksp->its = 0;
52: ksp->rnorm = dp;
53: PetscObjectGrantAccess(ksp);
54: KSPLogResidualHistory(ksp,dp);
55: KSPMonitor(ksp,0,dp);
56: (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
57: if (ksp->reason) return(0);
59: /* Make the initial Rp == R */
60: VecCopy(R,RP);
62: rhoold = 1.0;
63: alpha = 1.0;
64: omegaold = 1.0;
65: VecSet(P,0.0);
66: VecSet(V,0.0);
68: i=0;
69: do {
70: VecDot(R,RP,&rho); /* rho <- (r,rp) */
71: beta = (rho/rhoold) * (alpha/omegaold);
72: VecAXPBYPCZ(P,1.0,-omegaold*beta,beta,R,V); /* p <- r - omega * beta* v + beta * p */
73: KSP_PCApplyBAorAB(ksp,P,V,T); /* v <- K p */
74: VecDot(V,RP,&d1);
75: if (d1 == 0.0) SETERRQ(PETSC_ERR_PLIB,"Divide by zero");
76: alpha = rho / d1; /* a <- rho / (v,rp) */
77: VecWAXPY(S,-alpha,V,R); /* s <- r - a v */
78: KSP_PCApplyBAorAB(ksp,S,T,R);/* t <- K s */
79: VecDotNorm2(S,T,&d1,&d2);
80: if (d2 == 0.0) {
81: /* t is 0. if s is 0, then alpha v == r, and hence alpha p
82: may be our solution. Give it a try? */
83: VecDot(S,S,&d1);
84: if (d1 != 0.0) {
85: ksp->reason = KSP_DIVERGED_BREAKDOWN;
86: break;
87: }
88: VecAXPY(X,alpha,P); /* x <- x + a p */
89: PetscObjectTakeAccess(ksp);
90: ksp->its++;
91: ksp->rnorm = 0.0;
92: ksp->reason = KSP_CONVERGED_RTOL;
93: PetscObjectGrantAccess(ksp);
94: KSPLogResidualHistory(ksp,dp);
95: KSPMonitor(ksp,i+1,0.0);
96: break;
97: }
98: omega = d1 / d2; /* w <- (t's) / (t't) */
99: VecAXPBYPCZ(X,alpha,omega,1.0,P,S); /* x <- alpha * p + omega * s + x */
100: VecWAXPY(R,-omega,T,S); /* r <- s - w t */
101: if (ksp->normtype != KSP_NORM_NO && ksp->chknorm < i+2) {
102: VecNorm(R,NORM_2,&dp);
103: }
105: rhoold = rho;
106: omegaold = omega;
108: PetscObjectTakeAccess(ksp);
109: ksp->its++;
110: ksp->rnorm = dp;
111: PetscObjectGrantAccess(ksp);
112: KSPLogResidualHistory(ksp,dp);
113: KSPMonitor(ksp,i+1,dp);
114: (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
115: if (ksp->reason) break;
116: if (rho == 0.0) {
117: ksp->reason = KSP_DIVERGED_BREAKDOWN;
118: break;
119: }
120: i++;
121: } while (i<ksp->max_it);
123: if (i >= ksp->max_it) {
124: ksp->reason = KSP_DIVERGED_ITS;
125: }
127: KSPUnwindPreconditioner(ksp,X,T);
128: return(0);
129: }
131: /*MC
132: KSPBCGS - Implements the BiCGStab (Stabilized version of BiConjugate Gradient Squared) method.
134: Options Database Keys:
135: . see KSPSolve()
137: Level: beginner
139: Notes: Reference: van der Vorst, SIAM J. Sci. Stat. Comput., 1992.
140: See KSPBCGSL for additional stabilization
143: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPBICG, KSPBCGSL
144: M*/
148: PetscErrorCode KSPCreate_BCGS(KSP ksp)
149: {
151: ksp->data = (void*)0;
152: ksp->pc_side = PC_LEFT;
153: ksp->ops->setup = KSPSetUp_BCGS;
154: ksp->ops->solve = KSPSolve_BCGS;
155: ksp->ops->destroy = KSPDefaultDestroy;
156: ksp->ops->buildsolution = KSPDefaultBuildSolution;
157: ksp->ops->buildresidual = KSPDefaultBuildResidual;
158: ksp->ops->setfromoptions = 0;
159: ksp->ops->view = 0;
160: return(0);
161: }