xref: /petsc/src/ksp/pc/impls/pbjacobi/pbjacobi.c (revision 4877389942e948c907dd43b3a2d019cdca14ea44)
1dba47a55SKris Buschelman 
24b9ad928SBarry Smith /*
34b9ad928SBarry Smith    Include files needed for the PBJacobi preconditioner:
44b9ad928SBarry Smith      pcimpl.h - private include file intended for use by all preconditioners
54b9ad928SBarry Smith */
64b9ad928SBarry Smith 
7af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
84b9ad928SBarry Smith 
94b9ad928SBarry Smith /*
104b9ad928SBarry Smith    Private context (data structure) for the PBJacobi preconditioner.
114b9ad928SBarry Smith */
124b9ad928SBarry Smith typedef struct {
13713ccfa9SJed Brown   const MatScalar *diag;
14c1ac3661SBarry Smith   PetscInt         bs, mbs;
154b9ad928SBarry Smith } PC_PBJacobi;
164b9ad928SBarry Smith 
17d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_PBJacobi_1(PC pc, Vec x, Vec y)
18d71ae5a4SJacob Faibussowitsch {
19bbead8a2SBarry Smith   PC_PBJacobi       *jac = (PC_PBJacobi *)pc->data;
20bbead8a2SBarry Smith   PetscInt           i, m = jac->mbs;
21bbead8a2SBarry Smith   const MatScalar   *diag = jac->diag;
22bbead8a2SBarry Smith   const PetscScalar *xx;
23bbead8a2SBarry Smith   PetscScalar       *yy;
24bbead8a2SBarry Smith 
25bbead8a2SBarry Smith   PetscFunctionBegin;
269566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
279566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yy));
282fa5cd67SKarl Rupp   for (i = 0; i < m; i++) yy[i] = diag[i] * xx[i];
299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
309566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yy));
319566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(m));
32bbead8a2SBarry Smith   PetscFunctionReturn(0);
33bbead8a2SBarry Smith }
34bbead8a2SBarry Smith 
35d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_PBJacobi_2(PC pc, Vec x, Vec y)
36d71ae5a4SJacob Faibussowitsch {
374b9ad928SBarry Smith   PC_PBJacobi       *jac = (PC_PBJacobi *)pc->data;
38c1ac3661SBarry Smith   PetscInt           i, m = jac->mbs;
3985f4f44aSBarry Smith   const MatScalar   *diag = jac->diag;
400e666e78SToby Isaac   PetscScalar        x0, x1, *yy;
410e666e78SToby Isaac   const PetscScalar *xx;
424b9ad928SBarry Smith 
434b9ad928SBarry Smith   PetscFunctionBegin;
449566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
459566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yy));
464b9ad928SBarry Smith   for (i = 0; i < m; i++) {
479371c9d4SSatish Balay     x0            = xx[2 * i];
489371c9d4SSatish Balay     x1            = xx[2 * i + 1];
494b9ad928SBarry Smith     yy[2 * i]     = diag[0] * x0 + diag[2] * x1;
504b9ad928SBarry Smith     yy[2 * i + 1] = diag[1] * x0 + diag[3] * x1;
514b9ad928SBarry Smith     diag += 4;
524b9ad928SBarry Smith   }
539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yy));
559566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(6.0 * m));
564b9ad928SBarry Smith   PetscFunctionReturn(0);
574b9ad928SBarry Smith }
58d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_PBJacobi_3(PC pc, Vec x, Vec y)
59d71ae5a4SJacob Faibussowitsch {
604b9ad928SBarry Smith   PC_PBJacobi       *jac = (PC_PBJacobi *)pc->data;
61c1ac3661SBarry Smith   PetscInt           i, m = jac->mbs;
6285f4f44aSBarry Smith   const MatScalar   *diag = jac->diag;
630e666e78SToby Isaac   PetscScalar        x0, x1, x2, *yy;
640e666e78SToby Isaac   const PetscScalar *xx;
654b9ad928SBarry Smith 
664b9ad928SBarry Smith   PetscFunctionBegin;
679566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
689566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yy));
694b9ad928SBarry Smith   for (i = 0; i < m; i++) {
709371c9d4SSatish Balay     x0 = xx[3 * i];
719371c9d4SSatish Balay     x1 = xx[3 * i + 1];
729371c9d4SSatish Balay     x2 = xx[3 * i + 2];
732fa5cd67SKarl Rupp 
744b9ad928SBarry Smith     yy[3 * i]     = diag[0] * x0 + diag[3] * x1 + diag[6] * x2;
754b9ad928SBarry Smith     yy[3 * i + 1] = diag[1] * x0 + diag[4] * x1 + diag[7] * x2;
764b9ad928SBarry Smith     yy[3 * i + 2] = diag[2] * x0 + diag[5] * x1 + diag[8] * x2;
774b9ad928SBarry Smith     diag += 9;
784b9ad928SBarry Smith   }
799566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
809566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yy));
819566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(15.0 * m));
824b9ad928SBarry Smith   PetscFunctionReturn(0);
834b9ad928SBarry Smith }
84d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_PBJacobi_4(PC pc, Vec x, Vec y)
85d71ae5a4SJacob Faibussowitsch {
864b9ad928SBarry Smith   PC_PBJacobi       *jac = (PC_PBJacobi *)pc->data;
87c1ac3661SBarry Smith   PetscInt           i, m = jac->mbs;
8885f4f44aSBarry Smith   const MatScalar   *diag = jac->diag;
890e666e78SToby Isaac   PetscScalar        x0, x1, x2, x3, *yy;
900e666e78SToby Isaac   const PetscScalar *xx;
914b9ad928SBarry Smith 
924b9ad928SBarry Smith   PetscFunctionBegin;
939566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
949566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yy));
954b9ad928SBarry Smith   for (i = 0; i < m; i++) {
969371c9d4SSatish Balay     x0 = xx[4 * i];
979371c9d4SSatish Balay     x1 = xx[4 * i + 1];
989371c9d4SSatish Balay     x2 = xx[4 * i + 2];
999371c9d4SSatish Balay     x3 = xx[4 * i + 3];
1002fa5cd67SKarl Rupp 
1014b9ad928SBarry Smith     yy[4 * i]     = diag[0] * x0 + diag[4] * x1 + diag[8] * x2 + diag[12] * x3;
1024b9ad928SBarry Smith     yy[4 * i + 1] = diag[1] * x0 + diag[5] * x1 + diag[9] * x2 + diag[13] * x3;
1034b9ad928SBarry Smith     yy[4 * i + 2] = diag[2] * x0 + diag[6] * x1 + diag[10] * x2 + diag[14] * x3;
1044b9ad928SBarry Smith     yy[4 * i + 3] = diag[3] * x0 + diag[7] * x1 + diag[11] * x2 + diag[15] * x3;
1054b9ad928SBarry Smith     diag += 16;
1064b9ad928SBarry Smith   }
1079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
1089566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yy));
1099566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(28.0 * m));
1104b9ad928SBarry Smith   PetscFunctionReturn(0);
1114b9ad928SBarry Smith }
112d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_PBJacobi_5(PC pc, Vec x, Vec y)
113d71ae5a4SJacob Faibussowitsch {
1144b9ad928SBarry Smith   PC_PBJacobi       *jac = (PC_PBJacobi *)pc->data;
115c1ac3661SBarry Smith   PetscInt           i, m = jac->mbs;
11685f4f44aSBarry Smith   const MatScalar   *diag = jac->diag;
1170e666e78SToby Isaac   PetscScalar        x0, x1, x2, x3, x4, *yy;
1180e666e78SToby Isaac   const PetscScalar *xx;
1194b9ad928SBarry Smith 
1204b9ad928SBarry Smith   PetscFunctionBegin;
1219566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
1229566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yy));
1234b9ad928SBarry Smith   for (i = 0; i < m; i++) {
1249371c9d4SSatish Balay     x0 = xx[5 * i];
1259371c9d4SSatish Balay     x1 = xx[5 * i + 1];
1269371c9d4SSatish Balay     x2 = xx[5 * i + 2];
1279371c9d4SSatish Balay     x3 = xx[5 * i + 3];
1289371c9d4SSatish Balay     x4 = xx[5 * i + 4];
1292fa5cd67SKarl Rupp 
1304b9ad928SBarry Smith     yy[5 * i]     = diag[0] * x0 + diag[5] * x1 + diag[10] * x2 + diag[15] * x3 + diag[20] * x4;
1314b9ad928SBarry Smith     yy[5 * i + 1] = diag[1] * x0 + diag[6] * x1 + diag[11] * x2 + diag[16] * x3 + diag[21] * x4;
1324b9ad928SBarry Smith     yy[5 * i + 2] = diag[2] * x0 + diag[7] * x1 + diag[12] * x2 + diag[17] * x3 + diag[22] * x4;
1334b9ad928SBarry Smith     yy[5 * i + 3] = diag[3] * x0 + diag[8] * x1 + diag[13] * x2 + diag[18] * x3 + diag[23] * x4;
1344b9ad928SBarry Smith     yy[5 * i + 4] = diag[4] * x0 + diag[9] * x1 + diag[14] * x2 + diag[19] * x3 + diag[24] * x4;
1354b9ad928SBarry Smith     diag += 25;
1364b9ad928SBarry Smith   }
1379566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
1389566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yy));
1399566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(45.0 * m));
1404b9ad928SBarry Smith   PetscFunctionReturn(0);
1414b9ad928SBarry Smith }
142d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_PBJacobi_6(PC pc, Vec x, Vec y)
143d71ae5a4SJacob Faibussowitsch {
1440e1b4bd6SMark F. Adams   PC_PBJacobi       *jac = (PC_PBJacobi *)pc->data;
1450e1b4bd6SMark F. Adams   PetscInt           i, m = jac->mbs;
1460e1b4bd6SMark F. Adams   const MatScalar   *diag = jac->diag;
1470e666e78SToby Isaac   PetscScalar        x0, x1, x2, x3, x4, x5, *yy;
1480e666e78SToby Isaac   const PetscScalar *xx;
1490e1b4bd6SMark F. Adams 
1500e1b4bd6SMark F. Adams   PetscFunctionBegin;
1519566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
1529566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yy));
1530e1b4bd6SMark F. Adams   for (i = 0; i < m; i++) {
1549371c9d4SSatish Balay     x0 = xx[6 * i];
1559371c9d4SSatish Balay     x1 = xx[6 * i + 1];
1569371c9d4SSatish Balay     x2 = xx[6 * i + 2];
1579371c9d4SSatish Balay     x3 = xx[6 * i + 3];
1589371c9d4SSatish Balay     x4 = xx[6 * i + 4];
1599371c9d4SSatish Balay     x5 = xx[6 * i + 5];
1602fa5cd67SKarl Rupp 
1610e1b4bd6SMark F. Adams     yy[6 * i]     = diag[0] * x0 + diag[6] * x1 + diag[12] * x2 + diag[18] * x3 + diag[24] * x4 + diag[30] * x5;
1620e1b4bd6SMark F. Adams     yy[6 * i + 1] = diag[1] * x0 + diag[7] * x1 + diag[13] * x2 + diag[19] * x3 + diag[25] * x4 + diag[31] * x5;
1630e1b4bd6SMark F. Adams     yy[6 * i + 2] = diag[2] * x0 + diag[8] * x1 + diag[14] * x2 + diag[20] * x3 + diag[26] * x4 + diag[32] * x5;
1640e1b4bd6SMark F. Adams     yy[6 * i + 3] = diag[3] * x0 + diag[9] * x1 + diag[15] * x2 + diag[21] * x3 + diag[27] * x4 + diag[33] * x5;
1650e1b4bd6SMark F. Adams     yy[6 * i + 4] = diag[4] * x0 + diag[10] * x1 + diag[16] * x2 + diag[22] * x3 + diag[28] * x4 + diag[34] * x5;
1660e1b4bd6SMark F. Adams     yy[6 * i + 5] = diag[5] * x0 + diag[11] * x1 + diag[17] * x2 + diag[23] * x3 + diag[29] * x4 + diag[35] * x5;
1670e1b4bd6SMark F. Adams     diag += 36;
1680e1b4bd6SMark F. Adams   }
1699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
1709566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yy));
1719566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(66.0 * m));
1720e1b4bd6SMark F. Adams   PetscFunctionReturn(0);
1730e1b4bd6SMark F. Adams }
174d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_PBJacobi_7(PC pc, Vec x, Vec y)
175d71ae5a4SJacob Faibussowitsch {
1760a4ca348SMatthew G Knepley   PC_PBJacobi       *jac = (PC_PBJacobi *)pc->data;
1770a4ca348SMatthew G Knepley   PetscInt           i, m = jac->mbs;
1780a4ca348SMatthew G Knepley   const MatScalar   *diag = jac->diag;
1790e666e78SToby Isaac   PetscScalar        x0, x1, x2, x3, x4, x5, x6, *yy;
1800e666e78SToby Isaac   const PetscScalar *xx;
1810a4ca348SMatthew G Knepley 
1820a4ca348SMatthew G Knepley   PetscFunctionBegin;
1839566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
1849566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yy));
1850a4ca348SMatthew G Knepley   for (i = 0; i < m; i++) {
1869371c9d4SSatish Balay     x0 = xx[7 * i];
1879371c9d4SSatish Balay     x1 = xx[7 * i + 1];
1889371c9d4SSatish Balay     x2 = xx[7 * i + 2];
1899371c9d4SSatish Balay     x3 = xx[7 * i + 3];
1909371c9d4SSatish Balay     x4 = xx[7 * i + 4];
1919371c9d4SSatish Balay     x5 = xx[7 * i + 5];
1929371c9d4SSatish Balay     x6 = xx[7 * i + 6];
1932fa5cd67SKarl Rupp 
1940a4ca348SMatthew G Knepley     yy[7 * i]     = diag[0] * x0 + diag[7] * x1 + diag[14] * x2 + diag[21] * x3 + diag[28] * x4 + diag[35] * x5 + diag[42] * x6;
1950a4ca348SMatthew G Knepley     yy[7 * i + 1] = diag[1] * x0 + diag[8] * x1 + diag[15] * x2 + diag[22] * x3 + diag[29] * x4 + diag[36] * x5 + diag[43] * x6;
1960a4ca348SMatthew G Knepley     yy[7 * i + 2] = diag[2] * x0 + diag[9] * x1 + diag[16] * x2 + diag[23] * x3 + diag[30] * x4 + diag[37] * x5 + diag[44] * x6;
1970a4ca348SMatthew G Knepley     yy[7 * i + 3] = diag[3] * x0 + diag[10] * x1 + diag[17] * x2 + diag[24] * x3 + diag[31] * x4 + diag[38] * x5 + diag[45] * x6;
1980a4ca348SMatthew G Knepley     yy[7 * i + 4] = diag[4] * x0 + diag[11] * x1 + diag[18] * x2 + diag[25] * x3 + diag[32] * x4 + diag[39] * x5 + diag[46] * x6;
1990a4ca348SMatthew G Knepley     yy[7 * i + 5] = diag[5] * x0 + diag[12] * x1 + diag[19] * x2 + diag[26] * x3 + diag[33] * x4 + diag[40] * x5 + diag[47] * x6;
2000a4ca348SMatthew G Knepley     yy[7 * i + 6] = diag[6] * x0 + diag[13] * x1 + diag[20] * x2 + diag[27] * x3 + diag[34] * x4 + diag[41] * x5 + diag[48] * x6;
2010a4ca348SMatthew G Knepley     diag += 49;
2020a4ca348SMatthew G Knepley   }
2039566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
2049566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yy));
2059566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(91.0 * m)); /* 2*bs2 - bs */
2060a4ca348SMatthew G Knepley   PetscFunctionReturn(0);
2070a4ca348SMatthew G Knepley }
208d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_PBJacobi_N(PC pc, Vec x, Vec y)
209d71ae5a4SJacob Faibussowitsch {
210f79daf70SMark Lohry   PC_PBJacobi       *jac = (PC_PBJacobi *)pc->data;
211f79daf70SMark Lohry   PetscInt           i, ib, jb;
212f79daf70SMark Lohry   const PetscInt     m    = jac->mbs;
213f79daf70SMark Lohry   const PetscInt     bs   = jac->bs;
214f79daf70SMark Lohry   const MatScalar   *diag = jac->diag;
215f79daf70SMark Lohry   PetscScalar       *yy;
216f79daf70SMark Lohry   const PetscScalar *xx;
217f79daf70SMark Lohry 
218f79daf70SMark Lohry   PetscFunctionBegin;
2199566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
2209566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yy));
221f79daf70SMark Lohry   for (i = 0; i < m; i++) {
222f79daf70SMark Lohry     for (ib = 0; ib < bs; ib++) {
223f79daf70SMark Lohry       PetscScalar rowsum = 0;
224ad540459SPierre Jolivet       for (jb = 0; jb < bs; jb++) rowsum += diag[ib + jb * bs] * xx[bs * i + jb];
225f79daf70SMark Lohry       yy[bs * i + ib] = rowsum;
226f79daf70SMark Lohry     }
227f79daf70SMark Lohry     diag += bs * bs;
228f79daf70SMark Lohry   }
2299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
2309566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yy));
2319566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops((2.0 * bs * bs - bs) * m)); /* 2*bs2 - bs */
232f79daf70SMark Lohry   PetscFunctionReturn(0);
233f79daf70SMark Lohry }
234faff7d23SJed Brown 
235d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_PBJacobi_N(PC pc, Vec x, Vec y)
236d71ae5a4SJacob Faibussowitsch {
2378cd12025SJed Brown   PC_PBJacobi       *jac = (PC_PBJacobi *)pc->data;
2388cd12025SJed Brown   PetscInt           i, j, k, m = jac->mbs, bs = jac->bs;
2398cd12025SJed Brown   const MatScalar   *diag = jac->diag;
2408cd12025SJed Brown   const PetscScalar *xx;
2418cd12025SJed Brown   PetscScalar       *yy;
2428cd12025SJed Brown 
2438cd12025SJed Brown   PetscFunctionBegin;
2449566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
2459566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yy));
2468cd12025SJed Brown   for (i = 0; i < m; i++) {
2478cd12025SJed Brown     for (j = 0; j < bs; j++) yy[i * bs + j] = 0.;
2488cd12025SJed Brown     for (j = 0; j < bs; j++) {
249ad540459SPierre Jolivet       for (k = 0; k < bs; k++) yy[i * bs + k] += diag[k * bs + j] * xx[i * bs + j];
2508cd12025SJed Brown     }
2518cd12025SJed Brown     diag += bs * bs;
2528cd12025SJed Brown   }
2539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
2549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yy));
2559566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(m * bs * (2 * bs - 1)));
2568cd12025SJed Brown   PetscFunctionReturn(0);
2578cd12025SJed Brown }
2588cd12025SJed Brown 
259d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_PBJacobi(PC pc)
260d71ae5a4SJacob Faibussowitsch {
2614b9ad928SBarry Smith   PC_PBJacobi   *jac = (PC_PBJacobi *)pc->data;
2624b9ad928SBarry Smith   Mat            A   = pc->pmat;
26300e125f8SBarry Smith   MatFactorError err;
264539c167fSBarry Smith   PetscInt       nlocal;
2654b9ad928SBarry Smith 
2664b9ad928SBarry Smith   PetscFunctionBegin;
2679566063dSJacob Faibussowitsch   PetscCall(MatInvertBlockDiagonal(A, &jac->diag));
2689566063dSJacob Faibussowitsch   PetscCall(MatFactorGetError(A, &err));
26900e125f8SBarry Smith   if (err) pc->failedreason = (PCFailedReason)err;
270d7806518SHong Zhang 
2719566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A, &jac->bs));
2729566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &nlocal, NULL));
273539c167fSBarry Smith   jac->mbs = nlocal / jac->bs;
274521d7252SBarry Smith   switch (jac->bs) {
275d71ae5a4SJacob Faibussowitsch   case 1:
276d71ae5a4SJacob Faibussowitsch     pc->ops->apply = PCApply_PBJacobi_1;
277d71ae5a4SJacob Faibussowitsch     break;
278d71ae5a4SJacob Faibussowitsch   case 2:
279d71ae5a4SJacob Faibussowitsch     pc->ops->apply = PCApply_PBJacobi_2;
280d71ae5a4SJacob Faibussowitsch     break;
281d71ae5a4SJacob Faibussowitsch   case 3:
282d71ae5a4SJacob Faibussowitsch     pc->ops->apply = PCApply_PBJacobi_3;
283d71ae5a4SJacob Faibussowitsch     break;
284d71ae5a4SJacob Faibussowitsch   case 4:
285d71ae5a4SJacob Faibussowitsch     pc->ops->apply = PCApply_PBJacobi_4;
286d71ae5a4SJacob Faibussowitsch     break;
287d71ae5a4SJacob Faibussowitsch   case 5:
288d71ae5a4SJacob Faibussowitsch     pc->ops->apply = PCApply_PBJacobi_5;
289d71ae5a4SJacob Faibussowitsch     break;
290d71ae5a4SJacob Faibussowitsch   case 6:
291d71ae5a4SJacob Faibussowitsch     pc->ops->apply = PCApply_PBJacobi_6;
292d71ae5a4SJacob Faibussowitsch     break;
293d71ae5a4SJacob Faibussowitsch   case 7:
294d71ae5a4SJacob Faibussowitsch     pc->ops->apply = PCApply_PBJacobi_7;
295d71ae5a4SJacob Faibussowitsch     break;
296d71ae5a4SJacob Faibussowitsch   default:
297d71ae5a4SJacob Faibussowitsch     pc->ops->apply = PCApply_PBJacobi_N;
298d71ae5a4SJacob Faibussowitsch     break;
2994b9ad928SBarry Smith   }
3008cd12025SJed Brown   pc->ops->applytranspose = PCApplyTranspose_PBJacobi_N;
3014b9ad928SBarry Smith   PetscFunctionReturn(0);
3024b9ad928SBarry Smith }
303f1580f4eSBarry Smith 
304d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_PBJacobi(PC pc)
305d71ae5a4SJacob Faibussowitsch {
3064b9ad928SBarry Smith   PetscFunctionBegin;
3074b9ad928SBarry Smith   /*
3084b9ad928SBarry Smith       Free the private data structure that was hanging off the PC
3094b9ad928SBarry Smith   */
3109566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
3114b9ad928SBarry Smith   PetscFunctionReturn(0);
3124b9ad928SBarry Smith }
313fa939db7SJed Brown 
314d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_PBJacobi(PC pc, PetscViewer viewer)
315d71ae5a4SJacob Faibussowitsch {
316fa939db7SJed Brown   PC_PBJacobi *jac = (PC_PBJacobi *)pc->data;
317fa939db7SJed Brown   PetscBool    iascii;
318fa939db7SJed Brown 
319fa939db7SJed Brown   PetscFunctionBegin;
3209566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
32148a46eb9SPierre Jolivet   if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, "  point-block size %" PetscInt_FMT "\n", jac->bs));
322fa939db7SJed Brown   PetscFunctionReturn(0);
323fa939db7SJed Brown }
324fa939db7SJed Brown 
32537a17b4dSBarry Smith /*MC
326422a814eSBarry Smith      PCPBJACOBI - Point block Jacobi preconditioner
327422a814eSBarry Smith 
32895452b02SPatrick Sanan    Notes:
329f1580f4eSBarry Smith     See `PCJACOBI` for diagonal Jacobi, `PCVPBJACOBI` for variable-size point block, and `PCBJACOBI` for large size blocks
330422a814eSBarry Smith 
331f1580f4eSBarry Smith    This works for `MATAIJ` and `MATBAIJ` matrices and uses the blocksize provided to the matrix
332422a814eSBarry Smith 
333422a814eSBarry Smith    Uses dense LU factorization with partial pivoting to invert the blocks; if a zero pivot
334422a814eSBarry Smith    is detected a PETSc error is generated.
335422a814eSBarry Smith 
33695452b02SPatrick Sanan    Developer Notes:
337f1580f4eSBarry Smith      This should support the `PCSetErrorIfFailure()` flag set to `PETSC_TRUE` to allow
338422a814eSBarry Smith      the factorization to continue even after a zero pivot is found resulting in a Nan and hence
339*48773899SPierre Jolivet      terminating `KSP` with a `KSP_DIVERGED_NANORINF` allowing
340422a814eSBarry Smith      a nonlinear solver/ODE integrator to recover without stopping the program as currently happens.
341422a814eSBarry Smith 
34290dfcc32SBarry Smith      Perhaps should provide an option that allows generation of a valid preconditioner
343f1580f4eSBarry Smith      even if a block is singular as the `PCJACOBI` does.
34437a17b4dSBarry Smith 
34537a17b4dSBarry Smith    Level: beginner
34637a17b4dSBarry Smith 
347db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCJACOBI`, `PCVPBJACOBI`, `PCBJACOBI`
34837a17b4dSBarry Smith M*/
34937a17b4dSBarry Smith 
350d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_PBJacobi(PC pc)
351d71ae5a4SJacob Faibussowitsch {
3524b9ad928SBarry Smith   PC_PBJacobi *jac;
3534b9ad928SBarry Smith 
3544b9ad928SBarry Smith   PetscFunctionBegin;
3554b9ad928SBarry Smith   /*
3564b9ad928SBarry Smith      Creates the private data structure for this preconditioner and
3574b9ad928SBarry Smith      attach it to the PC object.
3584b9ad928SBarry Smith   */
3594dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&jac));
3604b9ad928SBarry Smith   pc->data = (void *)jac;
3614b9ad928SBarry Smith 
3624b9ad928SBarry Smith   /*
3634b9ad928SBarry Smith      Initialize the pointers to vectors to ZERO; these will be used to store
3644b9ad928SBarry Smith      diagonal entries of the matrix for fast preconditioner application.
3654b9ad928SBarry Smith   */
3660a545947SLisandro Dalcin   jac->diag = NULL;
3674b9ad928SBarry Smith 
3684b9ad928SBarry Smith   /*
3694b9ad928SBarry Smith       Set the pointers for the functions that are provided above.
3704b9ad928SBarry Smith       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
3714b9ad928SBarry Smith       are called, they will automatically call these functions.  Note we
3724b9ad928SBarry Smith       choose not to provide a couple of these functions since they are
3734b9ad928SBarry Smith       not needed.
3744b9ad928SBarry Smith   */
3750a545947SLisandro Dalcin   pc->ops->apply               = NULL; /*set depending on the block size */
3760a545947SLisandro Dalcin   pc->ops->applytranspose      = NULL;
3774b9ad928SBarry Smith   pc->ops->setup               = PCSetUp_PBJacobi;
3784b9ad928SBarry Smith   pc->ops->destroy             = PCDestroy_PBJacobi;
3790a545947SLisandro Dalcin   pc->ops->setfromoptions      = NULL;
380fa939db7SJed Brown   pc->ops->view                = PCView_PBJacobi;
3810a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
3820a545947SLisandro Dalcin   pc->ops->applysymmetricleft  = NULL;
3830a545947SLisandro Dalcin   pc->ops->applysymmetricright = NULL;
3844b9ad928SBarry Smith   PetscFunctionReturn(0);
3854b9ad928SBarry Smith }
386