xref: /petsc/src/ksp/pc/impls/vpbjacobi/vpbjacobi.c (revision 4877389942e948c907dd43b3a2d019cdca14ea44)
1 #include <../src/ksp/pc/impls/vpbjacobi/vpbjacobi.h>
2 
3 static PetscErrorCode PCApply_VPBJacobi(PC pc, Vec x, Vec y)
4 {
5   PC_VPBJacobi      *jac = (PC_VPBJacobi *)pc->data;
6   PetscInt           i, ncnt = 0;
7   const MatScalar   *diag = jac->diag;
8   PetscInt           ib, jb, bs;
9   const PetscScalar *xx;
10   PetscScalar       *yy, x0, x1, x2, x3, x4, x5, x6;
11   PetscInt           nblocks;
12   const PetscInt    *bsizes;
13 
14   PetscFunctionBegin;
15   PetscCall(MatGetVariableBlockSizes(pc->pmat, &nblocks, &bsizes));
16   PetscCall(VecGetArrayRead(x, &xx));
17   PetscCall(VecGetArray(y, &yy));
18   for (i = 0; i < nblocks; i++) {
19     bs = bsizes[i];
20     switch (bs) {
21     case 1:
22       yy[ncnt] = *diag * xx[ncnt];
23       break;
24     case 2:
25       x0           = xx[ncnt];
26       x1           = xx[ncnt + 1];
27       yy[ncnt]     = diag[0] * x0 + diag[2] * x1;
28       yy[ncnt + 1] = diag[1] * x0 + diag[3] * x1;
29       break;
30     case 3:
31       x0           = xx[ncnt];
32       x1           = xx[ncnt + 1];
33       x2           = xx[ncnt + 2];
34       yy[ncnt]     = diag[0] * x0 + diag[3] * x1 + diag[6] * x2;
35       yy[ncnt + 1] = diag[1] * x0 + diag[4] * x1 + diag[7] * x2;
36       yy[ncnt + 2] = diag[2] * x0 + diag[5] * x1 + diag[8] * x2;
37       break;
38     case 4:
39       x0           = xx[ncnt];
40       x1           = xx[ncnt + 1];
41       x2           = xx[ncnt + 2];
42       x3           = xx[ncnt + 3];
43       yy[ncnt]     = diag[0] * x0 + diag[4] * x1 + diag[8] * x2 + diag[12] * x3;
44       yy[ncnt + 1] = diag[1] * x0 + diag[5] * x1 + diag[9] * x2 + diag[13] * x3;
45       yy[ncnt + 2] = diag[2] * x0 + diag[6] * x1 + diag[10] * x2 + diag[14] * x3;
46       yy[ncnt + 3] = diag[3] * x0 + diag[7] * x1 + diag[11] * x2 + diag[15] * x3;
47       break;
48     case 5:
49       x0           = xx[ncnt];
50       x1           = xx[ncnt + 1];
51       x2           = xx[ncnt + 2];
52       x3           = xx[ncnt + 3];
53       x4           = xx[ncnt + 4];
54       yy[ncnt]     = diag[0] * x0 + diag[5] * x1 + diag[10] * x2 + diag[15] * x3 + diag[20] * x4;
55       yy[ncnt + 1] = diag[1] * x0 + diag[6] * x1 + diag[11] * x2 + diag[16] * x3 + diag[21] * x4;
56       yy[ncnt + 2] = diag[2] * x0 + diag[7] * x1 + diag[12] * x2 + diag[17] * x3 + diag[22] * x4;
57       yy[ncnt + 3] = diag[3] * x0 + diag[8] * x1 + diag[13] * x2 + diag[18] * x3 + diag[23] * x4;
58       yy[ncnt + 4] = diag[4] * x0 + diag[9] * x1 + diag[14] * x2 + diag[19] * x3 + diag[24] * x4;
59       break;
60     case 6:
61       x0           = xx[ncnt];
62       x1           = xx[ncnt + 1];
63       x2           = xx[ncnt + 2];
64       x3           = xx[ncnt + 3];
65       x4           = xx[ncnt + 4];
66       x5           = xx[ncnt + 5];
67       yy[ncnt]     = diag[0] * x0 + diag[6] * x1 + diag[12] * x2 + diag[18] * x3 + diag[24] * x4 + diag[30] * x5;
68       yy[ncnt + 1] = diag[1] * x0 + diag[7] * x1 + diag[13] * x2 + diag[19] * x3 + diag[25] * x4 + diag[31] * x5;
69       yy[ncnt + 2] = diag[2] * x0 + diag[8] * x1 + diag[14] * x2 + diag[20] * x3 + diag[26] * x4 + diag[32] * x5;
70       yy[ncnt + 3] = diag[3] * x0 + diag[9] * x1 + diag[15] * x2 + diag[21] * x3 + diag[27] * x4 + diag[33] * x5;
71       yy[ncnt + 4] = diag[4] * x0 + diag[10] * x1 + diag[16] * x2 + diag[22] * x3 + diag[28] * x4 + diag[34] * x5;
72       yy[ncnt + 5] = diag[5] * x0 + diag[11] * x1 + diag[17] * x2 + diag[23] * x3 + diag[29] * x4 + diag[35] * x5;
73       break;
74     case 7:
75       x0           = xx[ncnt];
76       x1           = xx[ncnt + 1];
77       x2           = xx[ncnt + 2];
78       x3           = xx[ncnt + 3];
79       x4           = xx[ncnt + 4];
80       x5           = xx[ncnt + 5];
81       x6           = xx[ncnt + 6];
82       yy[ncnt]     = diag[0] * x0 + diag[7] * x1 + diag[14] * x2 + diag[21] * x3 + diag[28] * x4 + diag[35] * x5 + diag[42] * x6;
83       yy[ncnt + 1] = diag[1] * x0 + diag[8] * x1 + diag[15] * x2 + diag[22] * x3 + diag[29] * x4 + diag[36] * x5 + diag[43] * x6;
84       yy[ncnt + 2] = diag[2] * x0 + diag[9] * x1 + diag[16] * x2 + diag[23] * x3 + diag[30] * x4 + diag[37] * x5 + diag[44] * x6;
85       yy[ncnt + 3] = diag[3] * x0 + diag[10] * x1 + diag[17] * x2 + diag[24] * x3 + diag[31] * x4 + diag[38] * x5 + diag[45] * x6;
86       yy[ncnt + 4] = diag[4] * x0 + diag[11] * x1 + diag[18] * x2 + diag[25] * x3 + diag[32] * x4 + diag[39] * x5 + diag[46] * x6;
87       yy[ncnt + 5] = diag[5] * x0 + diag[12] * x1 + diag[19] * x2 + diag[26] * x3 + diag[33] * x4 + diag[40] * x5 + diag[47] * x6;
88       yy[ncnt + 6] = diag[6] * x0 + diag[13] * x1 + diag[20] * x2 + diag[27] * x3 + diag[34] * x4 + diag[41] * x5 + diag[48] * x6;
89       break;
90     default:
91       for (ib = 0; ib < bs; ib++) {
92         PetscScalar rowsum = 0;
93         for (jb = 0; jb < bs; jb++) rowsum += diag[ib + jb * bs] * xx[ncnt + jb];
94         yy[ncnt + ib] = rowsum;
95       }
96     }
97     ncnt += bsizes[i];
98     diag += bsizes[i] * bsizes[i];
99   }
100   PetscCall(VecRestoreArrayRead(x, &xx));
101   PetscCall(VecRestoreArray(y, &yy));
102   PetscFunctionReturn(0);
103 }
104 
105 PETSC_INTERN PetscErrorCode PCSetUp_VPBJacobi_Host(PC pc)
106 {
107   PC_VPBJacobi   *jac = (PC_VPBJacobi *)pc->data;
108   Mat             A   = pc->pmat;
109   MatFactorError  err;
110   PetscInt        i, nsize = 0, nlocal;
111   PetscInt        nblocks;
112   const PetscInt *bsizes;
113 
114   PetscFunctionBegin;
115   PetscCall(MatGetVariableBlockSizes(pc->pmat, &nblocks, &bsizes));
116   PetscCall(MatGetLocalSize(pc->pmat, &nlocal, NULL));
117   PetscCheck(!nlocal || nblocks, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call MatSetVariableBlockSizes() before using PCVPBJACOBI");
118   if (!jac->diag) {
119     for (i = 0; i < nblocks; i++) nsize += bsizes[i] * bsizes[i];
120     PetscCall(PetscMalloc1(nsize, &jac->diag));
121   }
122   PetscCall(MatInvertVariableBlockDiagonal(A, nblocks, bsizes, jac->diag));
123   PetscCall(MatFactorGetError(A, &err));
124   if (err) pc->failedreason = (PCFailedReason)err;
125   pc->ops->apply = PCApply_VPBJacobi;
126   PetscFunctionReturn(0);
127 }
128 
129 static PetscErrorCode PCSetUp_VPBJacobi(PC pc)
130 {
131   PetscFunctionBegin;
132   /* In PCCreate_VPBJacobi() pmat might have not been set, so we wait to the last minute to do the dispatch */
133 #if defined(PETSC_HAVE_CUDA)
134   PetscBool isCuda;
135   PetscCall(PetscObjectTypeCompareAny((PetscObject)pc->pmat, &isCuda, MATSEQAIJCUSPARSE, MATMPIAIJCUSPARSE, ""));
136 #endif
137 #if defined(PETSC_HAVE_KOKKOS_KERNELS)
138   PetscBool isKok;
139   PetscCall(PetscObjectTypeCompareAny((PetscObject)pc->pmat, &isKok, MATSEQAIJKOKKOS, MATMPIAIJKOKKOS, ""));
140 #endif
141 
142 #if defined(PETSC_HAVE_CUDA)
143   if (isCuda) PetscCall(PCSetUp_VPBJacobi_CUDA(pc));
144   else
145 #endif
146 #if defined(PETSC_HAVE_KOKKOS_KERNELS)
147     if (isKok)
148     PetscCall(PCSetUp_VPBJacobi_Kokkos(pc));
149   else
150 #endif
151   {
152     PetscCall(PCSetUp_VPBJacobi_Host(pc));
153   }
154   PetscFunctionReturn(0);
155 }
156 
157 PETSC_INTERN PetscErrorCode PCDestroy_VPBJacobi(PC pc)
158 {
159   PC_VPBJacobi *jac = (PC_VPBJacobi *)pc->data;
160 
161   PetscFunctionBegin;
162   /*
163       Free the private data structure that was hanging off the PC
164   */
165   PetscCall(PetscFree(jac->diag));
166   PetscCall(PetscFree(pc->data));
167   PetscFunctionReturn(0);
168 }
169 
170 /*MC
171      PCVPBJACOBI - Variable size point block Jacobi preconditioner
172 
173    Level: beginner
174 
175    Notes:
176      See `PCJACOBI` for point Jacobi preconditioning, `PCPBJACOBI` for fixed point block size, and `PCBJACOBI` for large size blocks
177 
178      This works for `MATAIJ` matrices
179 
180      Uses dense LU factorization with partial pivoting to invert the blocks; if a zero pivot
181      is detected a PETSc error is generated.
182 
183      One must call `MatSetVariableBlockSizes()` to use this preconditioner
184 
185    Developer Notes:
186      This should support the `PCSetErrorIfFailure()` flag set to `PETSC_TRUE` to allow
187      the factorization to continue even after a zero pivot is found resulting in a Nan and hence
188      terminating `KSP` with a `KSP_DIVERGED_NANORINF` allowing
189      a nonlinear solver/ODE integrator to recover without stopping the program as currently happens.
190 
191      Perhaps should provide an option that allows generation of a valid preconditioner
192      even if a block is singular as the `PCJACOBI` does.
193 
194 .seealso: `MatSetVariableBlockSizes()`, `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCJACOBI`, `PCPBJACOBI`, `PCBJACOBI`
195 M*/
196 
197 PETSC_EXTERN PetscErrorCode PCCreate_VPBJacobi(PC pc)
198 {
199   PC_VPBJacobi *jac;
200 
201   PetscFunctionBegin;
202   /*
203      Creates the private data structure for this preconditioner and
204      attach it to the PC object.
205   */
206   PetscCall(PetscNew(&jac));
207   pc->data = (void *)jac;
208 
209   /*
210      Initialize the pointers to vectors to ZERO; these will be used to store
211      diagonal entries of the matrix for fast preconditioner application.
212   */
213   jac->diag = NULL;
214 
215   /*
216       Set the pointers for the functions that are provided above.
217       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
218       are called, they will automatically call these functions.  Note we
219       choose not to provide a couple of these functions since they are
220       not needed.
221   */
222   pc->ops->apply               = PCApply_VPBJacobi;
223   pc->ops->applytranspose      = NULL;
224   pc->ops->setup               = PCSetUp_VPBJacobi;
225   pc->ops->destroy             = PCDestroy_VPBJacobi;
226   pc->ops->setfromoptions      = NULL;
227   pc->ops->applyrichardson     = NULL;
228   pc->ops->applysymmetricleft  = NULL;
229   pc->ops->applysymmetricright = NULL;
230   PetscFunctionReturn(0);
231 }
232