Actual source code: bicg.c
1: /*$Id: bicg.c,v 1.26 2001/03/23 23:23:46 balay Exp $*/
3: /*
4: This code implements the BiCG (BiConjugate Gradient) method
6: Contibuted 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 pres,diagonalscale;
33: Scalar 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: pres = ksp->use_pres;
44: maxit = ksp->max_it;
45: X = ksp->vec_sol;
46: B = ksp->vec_rhs;
47: Rl = ksp->work[0];
48: Zl = ksp->work[1];
49: Pl = ksp->work[2];
50: Rr = ksp->work[3];
51: Zr = ksp->work[4];
52: Pr = ksp->work[5];
54: PCGetOperators(ksp->B,&Amat,&Pmat,&pflag);
56: if (!ksp->guess_zero) {
57: KSP_MatMult(ksp,Amat,X,Rr); /* r <- b - Ax */
58: VecAYPX(&mone,B,Rr);
59: } else {
60: VecCopy(B,Rr); /* r <- b (x is 0) */
61: }
62: VecCopy(Rr,Rl);
63: KSP_PCApply(ksp,ksp->B,Rr,Zr); /* z <- Br */
64: VecConjugate(Rl);
65: KSP_PCApplyTranspose(ksp,ksp->B,Rl,Zl);
66: VecConjugate(Rl);
67: VecConjugate(Zl);
68: if (pres) {
69: VecNorm(Zr,NORM_2,&dp); /* dp <- z'*z */
70: } else {
71: VecNorm(Rr,NORM_2,&dp); /* dp <- r'*r */
72: }
73: (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
74: if (ksp->reason) {*its = 0; return(0);}
75: KSPMonitor(ksp,0,dp);
76: PetscObjectTakeAccess(ksp);
77: ksp->its = 0;
78: ksp->rnorm = dp;
79: PetscObjectGrantAccess(ksp);
80: KSPLogResidualHistory(ksp,dp);
82: for (i=0; i<maxit; i++) {
83: VecDot(Zr,Rl,&beta); /* beta <- r'z */
84: if (!i) {
85: if (beta == 0.0) {
86: ksp->reason = KSP_DIVERGED_BREAKDOWN_BICG;
87: *its = 0;
88: return(0);
89: }
90: VecCopy(Zr,Pr); /* p <- z */
91: VecCopy(Zl,Pl);
92: } else {
93: b = beta/betaold;
94: VecAYPX(&b,Zr,Pr); /* p <- z + b* p */
95: b = PetscConj(b);
96: VecAYPX(&b,Zl,Pl);
97: }
98: betaold = beta;
99: KSP_MatMult(ksp,Amat,Pr,Zr); /* z <- Kp */
100: VecConjugate(Pl);
101: KSP_MatMultTranspose(ksp,Amat,Pl,Zl);
102: VecConjugate(Pl);
103: VecConjugate(Zl);
104: VecDot(Zr,Pl,&dpi); /* dpi <- z'p */
105: a = beta/dpi; /* a = beta/p'z */
106: VecAXPY(&a,Pr,X); /* x <- x + ap */
107: ma = -a;
108: VecAXPY(&ma,Zr,Rr);CHKERRQ(ierr)
109: ma = PetscConj(ma);
110: VecAXPY(&ma,Zl,Rl);
111: if (pres) {
112: KSP_PCApply(ksp,ksp->B,Rr,Zr); /* z <- Br */
113: VecConjugate(Rl);
114: KSP_PCApplyTranspose(ksp,ksp->B,Rl,Zl);
115: VecConjugate(Rl);
116: VecConjugate(Zl);
117: VecNorm(Zr,NORM_2,&dp); /* dp <- z'*z */
118: } else {
119: VecNorm(Rr,NORM_2,&dp); /* dp <- r'*r */
120: }
121: PetscObjectTakeAccess(ksp);
122: ksp->its = i+1;
123: ksp->rnorm = dp;
124: PetscObjectGrantAccess(ksp);
125: KSPLogResidualHistory(ksp,dp);
126: KSPMonitor(ksp,i+1,dp);
127: (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
128: if (ksp->reason) break;
129: if (!pres) {
130: KSP_PCApply(ksp,ksp->B,Rr,Zr); /* z <- Br */
131: VecConjugate(Rl);
132: KSP_PCApplyTranspose(ksp,ksp->B,Rl,Zl);
133: VecConjugate(Rl);
134: VecConjugate(Zl);
135: }
136: }
137: if (i == maxit) {i--; ksp->its--;ksp->reason = KSP_DIVERGED_ITS;}
138: *its = i+1;
139: return(0);
140: }
142: int KSPDestroy_BiCG(KSP ksp)
143: {
147: KSPDefaultFreeWork(ksp);
148: return(0);
149: }
151: EXTERN_C_BEGIN
152: int KSPCreate_BiCG(KSP ksp)
153: {
155: ksp->data = (void*)0;
156: ksp->pc_side = PC_LEFT;
157: ksp->calc_res = PETSC_TRUE;
158: ksp->ops->setup = KSPSetUp_BiCG;
159: ksp->ops->solve = KSPSolve_BiCG;
160: ksp->ops->destroy = KSPDestroy_BiCG;
161: ksp->ops->view = 0;
162: ksp->ops->setfromoptions = 0;
163: ksp->ops->buildsolution = KSPDefaultBuildSolution;
164: ksp->ops->buildresidual = KSPDefaultBuildResidual;
166: return(0);
167: }
168: EXTERN_C_END