Actual source code: bcgs.c
1: /*$Id: bcgs.c,v 1.77 2001/03/23 23:23:33 balay Exp $*/
3: /*
4: This code implements the BiCGStab (Stabilized version of BiConjugate
5: Gradient Squared) method. Reference: van der Vorst, SIAM J. Sci. Stat. Comput., 1992.
7: Note that for the complex numbers version, the VecDot() arguments
8: within the code MUST remain in the order given for correct computation
9: of inner products.
10: */
11: #include "src/sles/ksp/kspimpl.h"
13: static int KSPSetUp_BCGS(KSP ksp)
14: {
18: if (ksp->pc_side == PC_SYMMETRIC) {
19: SETERRQ(PETSC_ERR_SUP,"no symmetric preconditioning for KSPBCGS");
20: }
21: KSPDefaultGetWork(ksp,7);
22: return(0);
23: }
25: static int KSPSolve_BCGS(KSP ksp,int *its)
26: {
27: int i,maxit,ierr;
28: Scalar rho,rhoold,alpha,beta,omega,omegaold,d1,d2,zero = 0.0,tmp;
29: Vec X,B,V,P,R,RP,T,S,BINVF;
30: PetscReal dp = 0.0;
34: maxit = ksp->max_it;
35: X = ksp->vec_sol;
36: B = ksp->vec_rhs;
37: R = ksp->work[0];
38: RP = ksp->work[1];
39: V = ksp->work[2];
40: T = ksp->work[3];
41: S = ksp->work[4];
42: P = ksp->work[5];
43: BINVF = ksp->work[6];
45: /* Compute initial preconditioned residual */
46: KSPInitialResidual(ksp,X,V,T,R,BINVF,B);
48: /* Test for nothing to do */
49: if (!ksp->avoidnorms) {
50: VecNorm(R,NORM_2,&dp);
51: }
52: PetscObjectTakeAccess(ksp);
53: ksp->its = 0;
54: ksp->rnorm = dp;
55: PetscObjectGrantAccess(ksp);
56: KSPLogResidualHistory(ksp,dp);
57: KSPMonitor(ksp,0,dp);
58: (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
59: if (ksp->reason) {*its = 0; return(0);}
61: /* Make the initial Rp == R */
62: VecCopy(R,RP);
64: rhoold = 1.0;
65: alpha = 1.0;
66: omegaold = 1.0;
67: VecSet(&zero,P);
68: VecSet(&zero,V);
70: for (i=0; i<maxit; i++) {
72: VecDot(R,RP,&rho); /* rho <- (r,rp) */
73: if (rho == 0.0) SETERRQ(PETSC_ERR_KSP_BRKDWN,"Breakdown, rho = r . rp = 0");
74: beta = (rho/rhoold) * (alpha/omegaold);
75: tmp = -omegaold; VecAXPY(&tmp,V,P); /* p <- p - w v */
76: VecAYPX(&beta,R,P); /* p <- r + p beta */
77: KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,P,V,T); /* v <- K p */
78: VecDot(V,RP,&d1);
79: alpha = rho / d1; tmp = -alpha; /* a <- rho / (v,rp) */
80: VecWAXPY(&tmp,V,R,S); /* s <- r - a v */
81: KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,S,T,R);/* t <- K s */
82: VecDot(S,T,&d1);
83: VecDot(T,T,&d2);
84: if (d2 == 0.0) {
85: /* t is 0. if s is 0, then alpha v == r, and hence alpha p
86: may be our solution. Give it a try? */
87: VecDot(S,S,&d1);
88: if (d1 != 0.0) SETERRQ(PETSC_ERR_KSP_BRKDWN,"Breakdown, da = s . s != 0");
89: VecAXPY(&alpha,P,X); /* x <- x + a p */
90: PetscObjectTakeAccess(ksp);
91: ksp->its++;
92: ksp->rnorm = 0.0;
93: ksp->reason = KSP_CONVERGED_RTOL;
94: PetscObjectGrantAccess(ksp);
95: KSPLogResidualHistory(ksp,dp);
96: KSPMonitor(ksp,i+1,0.0);
97: break;
98: }
99: omega = d1 / d2; /* w <- (t's) / (t't) */
100: ierr = VecAXPY(&alpha,P,X); /* x <- x + a p */
101: ierr = VecAXPY(&omega,S,X); /* x <- x + w s */
102: tmp = -omega;
103: ierr = VecWAXPY(&tmp,T,S,R); /* r <- s - w t */
104: if (!ksp->avoidnorms) {
105: VecNorm(R,NORM_2,&dp);
106: }
108: rhoold = rho;
109: omegaold = omega;
111: PetscObjectTakeAccess(ksp);
112: ksp->its++;
113: ksp->rnorm = dp;
114: PetscObjectGrantAccess(ksp);
115: KSPLogResidualHistory(ksp,dp);
116: KSPMonitor(ksp,i+1,dp);
117: (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
118: if (ksp->reason) break;
119: }
120: if (i == maxit) {
121: ksp->reason = KSP_DIVERGED_ITS;
122: i--;
123: }
124: *its = i+1;
126: KSPUnwindPreconditioner(ksp,X,T);
127: return(0);
128: }
130: EXTERN_C_BEGIN
131: int KSPCreate_BCGS(KSP ksp)
132: {
134: ksp->data = (void*)0;
135: ksp->pc_side = PC_LEFT;
136: ksp->calc_res = PETSC_TRUE;
137: ksp->ops->setup = KSPSetUp_BCGS;
138: ksp->ops->solve = KSPSolve_BCGS;
139: ksp->ops->destroy = KSPDefaultDestroy;
140: ksp->ops->buildsolution = KSPDefaultBuildSolution;
141: ksp->ops->buildresidual = KSPDefaultBuildResidual;
142: ksp->ops->setfromoptions = 0;
143: ksp->ops->view = 0;
144: return(0);
145: }
146: EXTERN_C_END