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: #undef __FUNCT__
13: int KSPSetUp_BiCG(KSP ksp)
14: {
18: /* check user parameters and functions */
19: if (ksp->pc_side == PC_RIGHT) {
20: SETERRQ(2,"no right preconditioning for KSPBiCG");
21: } else if (ksp->pc_side == PC_SYMMETRIC) {
22: SETERRQ(2,"no symmetric preconditioning for KSPBiCG");
23: }
25: /* get work vectors from user code */
26: KSPDefaultGetWork(ksp,6);
28: return(0);
29: }
31: #undef __FUNCT__
33: int KSPSolve_BiCG(KSP ksp,int *its)
34: {
35: int ierr,i,maxit;
36: PetscTruth diagonalscale;
37: PetscScalar dpi,a=1.0,beta,betaold=1.0,b,mone=-1.0,ma;
38: PetscReal dp;
39: Vec X,B,Zl,Zr,Rl,Rr,Pl,Pr;
40: Mat Amat,Pmat;
41: MatStructure pflag;
44: ierr = PCDiagonalScale(ksp->B,&diagonalscale);
45: if (diagonalscale) SETERRQ1(1,"Krylov method %s does not support diagonal scaling",ksp->type_name);
47: maxit = ksp->max_it;
48: X = ksp->vec_sol;
49: B = ksp->vec_rhs;
50: Rl = ksp->work[0];
51: Zl = ksp->work[1];
52: Pl = ksp->work[2];
53: Rr = ksp->work[3];
54: Zr = ksp->work[4];
55: Pr = ksp->work[5];
57: PCGetOperators(ksp->B,&Amat,&Pmat,&pflag);
59: if (!ksp->guess_zero) {
60: KSP_MatMult(ksp,Amat,X,Rr); /* r <- b - Ax */
61: VecAYPX(&mone,B,Rr);
62: } else {
63: VecCopy(B,Rr); /* r <- b (x is 0) */
64: }
65: VecCopy(Rr,Rl);
66: KSP_PCApply(ksp,ksp->B,Rr,Zr); /* z <- Br */
67: VecConjugate(Rl);
68: KSP_PCApplyTranspose(ksp,ksp->B,Rl,Zl);
69: VecConjugate(Rl);
70: VecConjugate(Zl);
71: if (ksp->normtype == KSP_PRECONDITIONED_NORM) {
72: VecNorm(Zr,NORM_2,&dp); /* dp <- z'*z */
73: } else {
74: VecNorm(Rr,NORM_2,&dp); /* dp <- r'*r */
75: }
76: (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
77: if (ksp->reason) {*its = 0; return(0);}
78: KSPMonitor(ksp,0,dp);
79: PetscObjectTakeAccess(ksp);
80: ksp->its = 0;
81: ksp->rnorm = dp;
82: PetscObjectGrantAccess(ksp);
83: KSPLogResidualHistory(ksp,dp);
85: for (i=0; i<maxit; i++) {
86: VecDot(Zr,Rl,&beta); /* beta <- r'z */
87: if (!i) {
88: if (beta == 0.0) {
89: ksp->reason = KSP_DIVERGED_BREAKDOWN_BICG;
90: *its = 0;
91: return(0);
92: }
93: VecCopy(Zr,Pr); /* p <- z */
94: VecCopy(Zl,Pl);
95: } else {
96: b = beta/betaold;
97: VecAYPX(&b,Zr,Pr); /* p <- z + b* p */
98: b = PetscConj(b);
99: VecAYPX(&b,Zl,Pl);
100: }
101: betaold = beta;
102: KSP_MatMult(ksp,Amat,Pr,Zr); /* z <- Kp */
103: VecConjugate(Pl);
104: KSP_MatMultTranspose(ksp,Amat,Pl,Zl);
105: VecConjugate(Pl);
106: VecConjugate(Zl);
107: VecDot(Zr,Pl,&dpi); /* dpi <- z'p */
108: a = beta/dpi; /* a = beta/p'z */
109: VecAXPY(&a,Pr,X); /* x <- x + ap */
110: ma = -a;
111: VecAXPY(&ma,Zr,Rr);CHKERRQ(ierr)
112: ma = PetscConj(ma);
113: VecAXPY(&ma,Zl,Rl);
114: if (ksp->normtype == KSP_PRECONDITIONED_NORM) {
115: KSP_PCApply(ksp,ksp->B,Rr,Zr); /* z <- Br */
116: VecConjugate(Rl);
117: KSP_PCApplyTranspose(ksp,ksp->B,Rl,Zl);
118: VecConjugate(Rl);
119: VecConjugate(Zl);
120: VecNorm(Zr,NORM_2,&dp); /* dp <- z'*z */
121: } else {
122: VecNorm(Rr,NORM_2,&dp); /* dp <- r'*r */
123: }
124: PetscObjectTakeAccess(ksp);
125: ksp->its = i+1;
126: ksp->rnorm = dp;
127: PetscObjectGrantAccess(ksp);
128: KSPLogResidualHistory(ksp,dp);
129: KSPMonitor(ksp,i+1,dp);
130: (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
131: if (ksp->reason) break;
132: if (ksp->normtype == KSP_UNPRECONDITIONED_NORM) {
133: KSP_PCApply(ksp,ksp->B,Rr,Zr); /* z <- Br */
134: VecConjugate(Rl);
135: KSP_PCApplyTranspose(ksp,ksp->B,Rl,Zl);
136: VecConjugate(Rl);
137: VecConjugate(Zl);
138: }
139: }
140: if (i == maxit) {
141: ksp->reason = KSP_DIVERGED_ITS;
142: }
143: *its = ksp->its;
144: return(0);
145: }
147: #undef __FUNCT__
149: int KSPDestroy_BiCG(KSP ksp)
150: {
154: KSPDefaultFreeWork(ksp);
155: return(0);
156: }
158: EXTERN_C_BEGIN
159: #undef __FUNCT__
161: int KSPCreate_BiCG(KSP ksp)
162: {
164: ksp->data = (void*)0;
165: ksp->pc_side = PC_LEFT;
166: ksp->ops->setup = KSPSetUp_BiCG;
167: ksp->ops->solve = KSPSolve_BiCG;
168: ksp->ops->destroy = KSPDestroy_BiCG;
169: ksp->ops->view = 0;
170: ksp->ops->setfromoptions = 0;
171: ksp->ops->buildsolution = KSPDefaultBuildSolution;
172: ksp->ops->buildresidual = KSPDefaultBuildResidual;
174: return(0);
175: }
176: EXTERN_C_END