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