Actual source code: bcgs.c
1: /*$Id: bcgs.c,v 1.78 2001/08/07 03:03:49 balay Exp $*/
3: #include src/ksp/ksp/kspimpl.h
7: static int 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 int KSPSolve_BCGS(KSP ksp)
22: {
23: int i,ierr;
24: PetscScalar rho,rhoold,alpha,beta,omega,omegaold,d1,d2,zero = 0.0,tmp;
25: Vec X,B,V,P,R,RP,T,S;
26: PetscReal dp = 0.0;
30: X = ksp->vec_sol;
31: B = ksp->vec_rhs;
32: R = ksp->work[0];
33: RP = ksp->work[1];
34: V = ksp->work[2];
35: T = ksp->work[3];
36: S = ksp->work[4];
37: P = ksp->work[5];
39: /* Compute initial preconditioned residual */
40: KSPInitialResidual(ksp,X,V,T,R,B);
42: /* Test for nothing to do */
43: if (ksp->normtype != KSP_NO_NORM) {
44: VecNorm(R,NORM_2,&dp);
45: }
46: PetscObjectTakeAccess(ksp);
47: ksp->its = 0;
48: ksp->rnorm = dp;
49: PetscObjectGrantAccess(ksp);
50: KSPLogResidualHistory(ksp,dp);
51: KSPMonitor(ksp,0,dp);
52: (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
53: if (ksp->reason) return(0);
55: /* Make the initial Rp == R */
56: VecCopy(R,RP);
58: rhoold = 1.0;
59: alpha = 1.0;
60: omegaold = 1.0;
61: VecSet(&zero,P);
62: VecSet(&zero,V);
64: i=0;
65: do {
66: VecDot(R,RP,&rho); /* rho <- (r,rp) */
67: if (rho == 0.0) SETERRQ(PETSC_ERR_KSP_BRKDWN,"Breakdown, rho = r . rp = 0");
68: beta = (rho/rhoold) * (alpha/omegaold);
69: tmp = -omegaold; VecAXPY(&tmp,V,P); /* p <- p - w v */
70: VecAYPX(&beta,R,P); /* p <- r + p beta */
71: KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,P,V,T); /* v <- K p */
72: VecDot(V,RP,&d1);
73: alpha = rho / d1; tmp = -alpha; /* a <- rho / (v,rp) */
74: VecWAXPY(&tmp,V,R,S); /* s <- r - a v */
75: KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,S,T,R);/* t <- K s */
76: VecDot(S,T,&d1);
77: VecDot(T,T,&d2);
78: if (d2 == 0.0) {
79: /* t is 0. if s is 0, then alpha v == r, and hence alpha p
80: may be our solution. Give it a try? */
81: VecDot(S,S,&d1);
82: if (d1 != 0.0) SETERRQ(PETSC_ERR_KSP_BRKDWN,"Breakdown, da = s . s != 0");
83: VecAXPY(&alpha,P,X); /* x <- x + a p */
84: PetscObjectTakeAccess(ksp);
85: ksp->its++;
86: ksp->rnorm = 0.0;
87: ksp->reason = KSP_CONVERGED_RTOL;
88: PetscObjectGrantAccess(ksp);
89: KSPLogResidualHistory(ksp,dp);
90: KSPMonitor(ksp,i+1,0.0);
91: break;
92: }
93: omega = d1 / d2; /* w <- (t's) / (t't) */
94: VecAXPY(&alpha,P,X); /* x <- x + a p */
95: VecAXPY(&omega,S,X); /* x <- x + w s */
96: tmp = -omega;
97: VecWAXPY(&tmp,T,S,R); /* r <- s - w t */
98: if (ksp->normtype != KSP_NO_NORM) {
99: VecNorm(R,NORM_2,&dp);
100: }
102: rhoold = rho;
103: omegaold = omega;
105: PetscObjectTakeAccess(ksp);
106: ksp->its++;
107: ksp->rnorm = dp;
108: PetscObjectGrantAccess(ksp);
109: KSPLogResidualHistory(ksp,dp);
110: KSPMonitor(ksp,i+1,dp);
111: (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
112: if (ksp->reason) break;
113: i++;
114: } while (i<ksp->max_it);
116: if (i == ksp->max_it) {
117: ksp->reason = KSP_DIVERGED_ITS;
118: }
120: KSPUnwindPreconditioner(ksp,X,T);
121: return(0);
122: }
124: /*MC
125: KSPBCGS - Implements the BiCGStab (Stabilized version of BiConjugate
126: Gradient Squared) method. Reference: van der Vorst, SIAM J. Sci. Stat. Comput., 1992.
129: Options Database Keys:
130: . see KSPSolve()
132: Level: beginner
134: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPBICG
136: M*/
138: EXTERN_C_BEGIN
141: int KSPCreate_BCGS(KSP ksp)
142: {
144: ksp->data = (void*)0;
145: ksp->pc_side = PC_LEFT;
146: ksp->ops->setup = KSPSetUp_BCGS;
147: ksp->ops->solve = KSPSolve_BCGS;
148: ksp->ops->destroy = KSPDefaultDestroy;
149: ksp->ops->buildsolution = KSPDefaultBuildSolution;
150: ksp->ops->buildresidual = KSPDefaultBuildResidual;
151: ksp->ops->setfromoptions = 0;
152: ksp->ops->view = 0;
153: return(0);
154: }
155: EXTERN_C_END