Actual source code: ibcgs.c
2: #include <private/kspimpl.h>
6: static PetscErrorCode KSPSetUp_IBCGS(KSP ksp)
7: {
9: PetscBool diagonalscale;
12: PCGetDiagonalScale(ksp->pc,&diagonalscale);
13: if (diagonalscale) SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
14: KSPDefaultGetWork(ksp,9);
15: return(0);
16: }
18: /*
19: The code below "cheats" from PETSc style
20: 1) VecRestoreArray() is called immediately after VecGetArray() and the array values are still accessed; the reason for the immediate
21: restore is that Vec operations are done on some of the vectors during the solve and if we did not restore immediately it would
22: generate two VecGetArray() (the second one inside the Vec operation) calls without a restore between them.
23: 2) The vector operations on done directly on the arrays instead of with VecXXXX() calls
25: For clarity in the code we name single VECTORS with two names, for example, Rn_1 and R, but they actually always
26: the exact same memory. We do this with macro defines so that compiler won't think they are
27: two different variables.
29: */
30: #define Xn_1 Xn
31: #define xn_1 xn
32: #define Rn_1 Rn
33: #define rn_1 rn
34: #define Un_1 Un
35: #define un_1 un
36: #define Vn_1 Vn
37: #define vn_1 vn
38: #define Qn_1 Qn
39: #define qn_1 qn
40: #define Zn_1 Zn
41: #define zn_1 zn
44: static PetscErrorCode KSPSolve_IBCGS(KSP ksp)
45: {
47: PetscInt i,N;
48: PetscReal rnorm,rnormin = 0.0;
49: #if defined(PETSC_HAVE_MPI_LONG_DOUBLE) && !defined(PETSC_USE_COMPLEX)
50: /* Because of possible instabilities in the algorithm (as indicated by different residual histories for the same problem
51: on the same number of processes with different runs) we support computing the inner products using Intel's 80 bit arithematic
52: rather than just 64 bit. Thus we copy our double precision values into long doubles (hoping this keeps the 16 extra bits)
53: and tell MPI to do its ALlreduces with MPI_LONG_DOUBLE.
55: Note for developers that does not effect the code. Intel's long double is implemented by storing the 80 bits of extended double
56: precision into a 16 byte space (the rest of the space is ignored) */
57: long double insums[7],outsums[7];
58: #else
59: PetscScalar insums[7],outsums[7];
60: #endif
61: PetscScalar sigman_2, sigman_1, sigman, pin_1, pin, phin_1, phin,tmp1,tmp2;
62: PetscScalar taun_1, taun, rhon, alphan_1, alphan, omegan_1, omegan;
63: PetscScalar *PETSC_RESTRICT r0, *PETSC_RESTRICT rn, *PETSC_RESTRICT xn, *PETSC_RESTRICT f0, *PETSC_RESTRICT vn, *PETSC_RESTRICT zn, *PETSC_RESTRICT qn;
64: PetscScalar *PETSC_RESTRICT b, *PETSC_RESTRICT un;
65: /* the rest do not have to keep n_1 values */
66: PetscScalar kappan, thetan, etan, gamman, betan, deltan;
67: PetscScalar *PETSC_RESTRICT tn, *PETSC_RESTRICT sn;
68: Vec R0,Rn,Xn,F0,Vn,Zn,Qn,Tn,Sn,B,Un;
69: Mat A;
72: if (!ksp->vec_rhs->petscnative) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"Only coded for PETSc vectors");
74: PCGetOperators(ksp->pc,&A,PETSC_NULL,PETSC_NULL);
75: VecGetLocalSize(ksp->vec_sol,&N);
76: Xn = ksp->vec_sol;VecGetArray(Xn_1,(PetscScalar**)&xn_1);VecRestoreArray(Xn_1,(PetscScalar**)&xn_1);
77: B = ksp->vec_rhs;VecGetArray(B,(PetscScalar**)&b);VecRestoreArray(B,(PetscScalar**)&b);
78: R0 = ksp->work[0];VecGetArray(R0,(PetscScalar**)&r0);VecRestoreArray(R0,(PetscScalar**)&r0);
79: Rn = ksp->work[1];VecGetArray(Rn_1,(PetscScalar**)&rn_1);VecRestoreArray(Rn_1,(PetscScalar**)&rn_1);
80: Un = ksp->work[2];VecGetArray(Un_1,(PetscScalar**)&un_1);VecRestoreArray(Un_1,(PetscScalar**)&un_1);
81: F0 = ksp->work[3];VecGetArray(F0,(PetscScalar**)&f0);VecRestoreArray(F0,(PetscScalar**)&f0);
82: Vn = ksp->work[4];VecGetArray(Vn_1,(PetscScalar**)&vn_1);VecRestoreArray(Vn_1,(PetscScalar**)&vn_1);
83: Zn = ksp->work[5];VecGetArray(Zn_1,(PetscScalar**)&zn_1);VecRestoreArray(Zn_1,(PetscScalar**)&zn_1);
84: Qn = ksp->work[6];VecGetArray(Qn_1,(PetscScalar**)&qn_1);VecRestoreArray(Qn_1,(PetscScalar**)&qn_1);
85: Tn = ksp->work[7];VecGetArray(Tn,(PetscScalar**)&tn);VecRestoreArray(Tn,(PetscScalar**)&tn);
86: Sn = ksp->work[8];VecGetArray(Sn,(PetscScalar**)&sn);VecRestoreArray(Sn,(PetscScalar**)&sn);
88: /* r0 = rn_1 = b - A*xn_1; */
89: /* KSP_PCApplyBAorAB(ksp,Xn_1,Rn_1,Tn);
90: VecAYPX(Rn_1,-1.0,B); */
91: KSPInitialResidual(ksp,Xn_1,Tn,Sn,Rn_1,B);
93: VecNorm(Rn_1,NORM_2,&rnorm);
94: KSPMonitor(ksp,0,rnorm);
95: (*ksp->converged)(ksp,0,rnorm,&ksp->reason,ksp->cnvP);
96: if (ksp->reason) return(0);
98: VecCopy(Rn_1,R0);
100: /* un_1 = A*rn_1; */
101: KSP_PCApplyBAorAB(ksp,Rn_1,Un_1,Tn);
102:
103: /* f0 = A'*rn_1; */
104: if (ksp->pc_side == PC_RIGHT) { /* B' A' */
105: MatMultTranspose(A,R0,Tn);
106: PCApplyTranspose(ksp->pc,Tn,F0);
107: } else if (ksp->pc_side == PC_LEFT) { /* A' B' */
108: PCApplyTranspose(ksp->pc,R0,Tn);
109: MatMultTranspose(A,Tn,F0);
110: }
112: /*qn_1 = vn_1 = zn_1 = 0.0; */
113: VecSet(Qn_1,0.0);
114: VecSet(Vn_1,0.0);
115: VecSet(Zn_1,0.0);
117: sigman_2 = pin_1 = taun_1 = 0.0;
119: /* the paper says phin_1 should be initialized to zero, it is actually R0'R0 */
120: VecDot(R0,R0,&phin_1);
122: /* sigman_1 = rn_1'un_1 */
123: VecDot(R0,Un_1,&sigman_1);
125: alphan_1 = omegan_1 = 1.0;
127: for (ksp->its = 1; ksp->its<ksp->max_it+1; ksp->its++) {
128: rhon = phin_1 - omegan_1*sigman_2 + omegan_1*alphan_1*pin_1;
129: /* if (rhon == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED,"rhon is zero, iteration %D",n); */
130: if (ksp->its == 1) deltan = rhon;
131: else deltan = rhon/taun_1;
132: betan = deltan/omegan_1;
133: taun = sigman_1 + betan*taun_1 - deltan*pin_1;
134: if (taun == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED,"taun is zero, iteration %D",ksp->its);
135: alphan = rhon/taun;
136: PetscLogFlops(15.0);
138: /*
139: zn = alphan*rn_1 + (alphan/alphan_1)betan*zn_1 - alphan*deltan*vn_1
140: vn = un_1 + betan*vn_1 - deltan*qn_1
141: sn = rn_1 - alphan*vn
143: The algorithm in the paper is missing the alphan/alphan_1 term in the zn update
144: */
145: PetscLogEventBegin(VEC_Ops,0,0,0,0);
146: tmp1 = (alphan/alphan_1)*betan;
147: tmp2 = alphan*deltan;
148: for (i=0; i<N; i++) {
149: zn[i] = alphan*rn_1[i] + tmp1*zn_1[i] - tmp2*vn_1[i];
150: vn[i] = un_1[i] + betan*vn_1[i] - deltan*qn_1[i];
151: sn[i] = rn_1[i] - alphan*vn[i];
152: }
153: PetscLogFlops(3.0+11.0*N);
154: PetscLogEventEnd(VEC_Ops,0,0,0,0);
156: /*
157: qn = A*vn
158: */
159: KSP_PCApplyBAorAB(ksp,Vn,Qn,Tn);
161: /*
162: tn = un_1 - alphan*qn
163: */
164: VecWAXPY(Tn,-alphan,Qn,Un_1);
165:
167: /*
168: phin = r0'sn
169: pin = r0'qn
170: gamman = f0'sn
171: etan = f0'tn
172: thetan = sn'tn
173: kappan = tn'tn
174: */
175: PetscLogEventBegin(VEC_ReduceArithmetic,0,0,0,0);
176: phin = pin = gamman = etan = thetan = kappan = 0.0;
177: for (i=0; i<N; i++) {
178: phin += r0[i]*sn[i];
179: pin += r0[i]*qn[i];
180: gamman += f0[i]*sn[i];
181: etan += f0[i]*tn[i];
182: thetan += sn[i]*tn[i];
183: kappan += tn[i]*tn[i];
184: }
185: PetscLogFlops(12.0*N);
186: PetscLogEventEnd(VEC_ReduceArithmetic,0,0,0,0);
187: insums[0] = phin;
188: insums[1] = pin;
189: insums[2] = gamman;
190: insums[3] = etan;
191: insums[4] = thetan;
192: insums[5] = kappan;
193: insums[6] = rnormin;
194: PetscLogEventBarrierBegin(VEC_ReduceBarrier,0,0,0,0,((PetscObject)ksp)->comm);
195: #if defined(PETSC_HAVE_MPI_LONG_DOUBLE) && !defined(PETSC_USE_COMPLEX)
196: if (ksp->lagnorm && ksp->its > 1) {
197: MPI_Allreduce(insums,outsums,7,MPI_LONG_DOUBLE,MPI_SUM,((PetscObject)ksp)->comm);
198: } else {
199: MPI_Allreduce(insums,outsums,6,MPI_LONG_DOUBLE,MPI_SUM,((PetscObject)ksp)->comm);
200: }
201: #else
202: if (ksp->lagnorm && ksp->its > 1) {
203: MPI_Allreduce(insums,outsums,7,MPIU_SCALAR,MPI_SUM,((PetscObject)ksp)->comm);
204: } else {
205: MPI_Allreduce(insums,outsums,6,MPIU_SCALAR,MPI_SUM,((PetscObject)ksp)->comm);
206: }
207: #endif
208: PetscLogEventBarrierEnd(VEC_ReduceBarrier,0,0,0,0,((PetscObject)ksp)->comm);
209: phin = outsums[0];
210: pin = outsums[1];
211: gamman = outsums[2];
212: etan = outsums[3];
213: thetan = outsums[4];
214: kappan = outsums[5];
215: if (ksp->lagnorm && ksp->its > 1) rnorm = PetscSqrtReal(PetscRealPart(outsums[6]));
217: if (kappan == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED,"kappan is zero, iteration %D",ksp->its);
218: if (thetan == 0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED,"thetan is zero, iteration %D",ksp->its);
219: omegan = thetan/kappan;
220: sigman = gamman - omegan*etan;
222: /*
223: rn = sn - omegan*tn
224: xn = xn_1 + zn + omegan*sn
225: */
226: PetscLogEventBegin(VEC_Ops,0,0,0,0);
227: rnormin = 0.0;
228: for (i=0; i<N; i++) {
229: rn[i] = sn[i] - omegan*tn[i];
230: rnormin += PetscRealPart(PetscConj(rn[i])*rn[i]);
231: xn[i] += zn[i] + omegan*sn[i];
232: }
233: PetscLogFlops(7.0*N);
234: PetscLogEventEnd(VEC_Ops,0,0,0,0);
236: if (!ksp->lagnorm && ksp->chknorm < ksp->its) {
237: PetscLogEventBarrierBegin(VEC_ReduceBarrier,0,0,0,0,((PetscObject)ksp)->comm);
238: MPI_Allreduce(&rnormin,&rnorm,1,MPIU_REAL,MPIU_SUM,((PetscObject)ksp)->comm);
239: PetscLogEventBarrierEnd(VEC_ReduceBarrier,0,0,0,0,((PetscObject)ksp)->comm);
240: rnorm = PetscSqrtReal(rnorm);
241: }
243: /* Test for convergence */
244: KSPMonitor(ksp,ksp->its,rnorm);
245: (*ksp->converged)(ksp,ksp->its,rnorm,&ksp->reason,ksp->cnvP);
246: if (ksp->reason) break;
247:
248: /* un = A*rn */
249: KSP_PCApplyBAorAB(ksp,Rn,Un,Tn);
251: /* Update n-1 locations with n locations */
252: sigman_2 = sigman_1;
253: sigman_1 = sigman;
254: pin_1 = pin;
255: phin_1 = phin;
256: alphan_1 = alphan;
257: taun_1 = taun;
258: omegan_1 = omegan;
259: }
260: if (ksp->its >= ksp->max_it) {
261: ksp->reason = KSP_DIVERGED_ITS;
262: }
263: KSPUnwindPreconditioner(ksp,Xn,Tn);
264: return(0);
265: }
268: /*MC
269: KSPIBCGS - Implements the IBiCGStab (Improved Stabilized version of BiConjugate Gradient Squared) method
270: in an alternative form to have only a single global reduction operation instead of the usual 3 (or 4)
272: Options Database Keys:
273: . see KSPSolve()
275: Level: beginner
277: Notes: Supports left and right preconditioning
279: See KSPBCGSL for additional stabilization
281: Unlike the Bi-CG-stab algorithm, this requires one multiplication be the transpose of the operator
282: before the iteration starts.
284: The paper has two errors in the algorithm presented, they are fixed in the code in KSPSolve_IBCGS()
286: For maximum reduction in the number of global reduction operations, this solver should be used with
287: KSPSetLagNorm().
289: This is not supported for complex numbers.
291: Reference: The Improved BiCGStab Method for Large and Sparse Unsymmetric Linear Systems on Parallel Distributed Memory
292: Architectures. L. T. Yand and R. Brent, Proceedings of the Fifth International Conference on Algorithms and
293: Architectures for Parallel Processing, 2002, IEEE.
295: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPBICG, KSPBCGSL, KSPIBCGS, KSPSetLagNorm()
296: M*/
300: PetscErrorCode KSPCreate_IBCGS(KSP ksp)
301: {
305: ksp->data = (void*)0;
306: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
307: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,1);
308: ksp->ops->setup = KSPSetUp_IBCGS;
309: ksp->ops->solve = KSPSolve_IBCGS;
310: ksp->ops->destroy = KSPDefaultDestroy;
311: ksp->ops->buildsolution = KSPDefaultBuildSolution;
312: ksp->ops->buildresidual = KSPDefaultBuildResidual;
313: ksp->ops->setfromoptions = 0;
314: ksp->ops->view = 0;
315: #if defined(PETSC_USE_COMPLEX)
316: SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"This is not supported for complex numbers");
317: #endif
318: return(0);
319: }