Actual source code: bicg.c
1: /*$Id: bicg.c,v 1.28 2001/08/07 03:03:55 balay Exp $*/
3: /*
4: This code implements the BiCG (BiConjugate Gradient) method
6: Contributed by: Victor Eijkhout
8: */
9: #include src/sles/ksp/kspimpl.h
11: int KSPSetUp_BiCG(KSP ksp)
12: {
16: /* check user parameters and functions */
17: if (ksp->pc_side == PC_RIGHT) {
18: SETERRQ(2,"no right preconditioning for KSPBiCG");
19: } else if (ksp->pc_side == PC_SYMMETRIC) {
20: SETERRQ(2,"no symmetric preconditioning for KSPBiCG");
21: }
23: /* get work vectors from user code */
24: KSPDefaultGetWork(ksp,6);
26: return(0);
27: }
29: int KSPSolve_BiCG(KSP ksp,int *its)
30: {
31: int ierr,i,maxit;
32: PetscTruth diagonalscale;
33: PetscScalar dpi,a=1.0,beta,betaold=1.0,b,mone=-1.0,ma;
34: PetscReal dp;
35: Vec X,B,Zl,Zr,Rl,Rr,Pl,Pr;
36: Mat Amat,Pmat;
37: MatStructure pflag;
40: ierr = PCDiagonalScale(ksp->B,&diagonalscale);
41: if (diagonalscale) SETERRQ1(1,"Krylov method %s does not support diagonal scaling",ksp->type_name);
43: maxit = ksp->max_it;
44: X = ksp->vec_sol;
45: B = ksp->vec_rhs;
46: Rl = ksp->work[0];
47: Zl = ksp->work[1];
48: Pl = ksp->work[2];
49: Rr = ksp->work[3];
50: Zr = ksp->work[4];
51: Pr = ksp->work[5];
53: PCGetOperators(ksp->B,&Amat,&Pmat,&pflag);
55: if (!ksp->guess_zero) {
56: KSP_MatMult(ksp,Amat,X,Rr); /* r <- b - Ax */
57: VecAYPX(&mone,B,Rr);
58: } else {
59: VecCopy(B,Rr); /* r <- b (x is 0) */
60: }
61: VecCopy(Rr,Rl);
62: KSP_PCApply(ksp,ksp->B,Rr,Zr); /* z <- Br */
63: VecConjugate(Rl);
64: KSP_PCApplyTranspose(ksp,ksp->B,Rl,Zl);
65: VecConjugate(Rl);
66: VecConjugate(Zl);
67: if (ksp->normtype == KSP_PRECONDITIONED_NORM) {
68: VecNorm(Zr,NORM_2,&dp); /* dp <- z'*z */
69: } else {
70: VecNorm(Rr,NORM_2,&dp); /* dp <- r'*r */
71: }
72: (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
73: if (ksp->reason) {*its = 0; return(0);}
74: KSPMonitor(ksp,0,dp);
75: PetscObjectTakeAccess(ksp);
76: ksp->its = 0;
77: ksp->rnorm = dp;
78: PetscObjectGrantAccess(ksp);
79: KSPLogResidualHistory(ksp,dp);
81: for (i=0; i<maxit; i++) {
82: VecDot(Zr,Rl,&beta); /* beta <- r'z */
83: if (!i) {
84: if (beta == 0.0) {
85: ksp->reason = KSP_DIVERGED_BREAKDOWN_BICG;
86: *its = 0;
87: return(0);
88: }
89: VecCopy(Zr,Pr); /* p <- z */
90: VecCopy(Zl,Pl);
91: } else {
92: b = beta/betaold;
93: VecAYPX(&b,Zr,Pr); /* p <- z + b* p */
94: b = PetscConj(b);
95: VecAYPX(&b,Zl,Pl);
96: }
97: betaold = beta;
98: KSP_MatMult(ksp,Amat,Pr,Zr); /* z <- Kp */
99: VecConjugate(Pl);
100: KSP_MatMultTranspose(ksp,Amat,Pl,Zl);
101: VecConjugate(Pl);
102: VecConjugate(Zl);
103: VecDot(Zr,Pl,&dpi); /* dpi <- z'p */
104: a = beta/dpi; /* a = beta/p'z */
105: VecAXPY(&a,Pr,X); /* x <- x + ap */
106: ma = -a;
107: VecAXPY(&ma,Zr,Rr);CHKERRQ(ierr)
108: ma = PetscConj(ma);
109: VecAXPY(&ma,Zl,Rl);
110: if (ksp->normtype == KSP_PRECONDITIONED_NORM) {
111: KSP_PCApply(ksp,ksp->B,Rr,Zr); /* z <- Br */
112: VecConjugate(Rl);
113: KSP_PCApplyTranspose(ksp,ksp->B,Rl,Zl);
114: VecConjugate(Rl);
115: VecConjugate(Zl);
116: VecNorm(Zr,NORM_2,&dp); /* dp <- z'*z */
117: } else {
118: VecNorm(Rr,NORM_2,&dp); /* dp <- r'*r */
119: }
120: PetscObjectTakeAccess(ksp);
121: ksp->its = i+1;
122: ksp->rnorm = dp;
123: PetscObjectGrantAccess(ksp);
124: KSPLogResidualHistory(ksp,dp);
125: KSPMonitor(ksp,i+1,dp);
126: (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
127: if (ksp->reason) break;
128: if (ksp->normtype == KSP_UNPRECONDITIONED_NORM) {
129: KSP_PCApply(ksp,ksp->B,Rr,Zr); /* z <- Br */
130: VecConjugate(Rl);
131: KSP_PCApplyTranspose(ksp,ksp->B,Rl,Zl);
132: VecConjugate(Rl);
133: VecConjugate(Zl);
134: }
135: }
136: if (i == maxit) {i--; ksp->its--;ksp->reason = KSP_DIVERGED_ITS;}
137: *its = i+1;
138: return(0);
139: }
141: int KSPDestroy_BiCG(KSP ksp)
142: {
146: KSPDefaultFreeWork(ksp);
147: return(0);
148: }
150: EXTERN_C_BEGIN
151: int KSPCreate_BiCG(KSP ksp)
152: {
154: ksp->data = (void*)0;
155: ksp->pc_side = PC_LEFT;
156: ksp->ops->setup = KSPSetUp_BiCG;
157: ksp->ops->solve = KSPSolve_BiCG;
158: ksp->ops->destroy = KSPDestroy_BiCG;
159: ksp->ops->view = 0;
160: ksp->ops->setfromoptions = 0;
161: ksp->ops->buildsolution = KSPDefaultBuildSolution;
162: ksp->ops->buildresidual = KSPDefaultBuildResidual;
164: return(0);
165: }
166: EXTERN_C_END