xref: /petsc/src/ksp/pc/impls/bjacobi/bjacobi.c (revision e0f629dd96229b85b8d4fbf4a5ed4f197d232daa)
1 
2 /*
3    Defines a block Jacobi preconditioner.
4 */
5 
6 #include <../src/ksp/pc/impls/bjacobi/bjacobi.h> /*I "petscpc.h" I*/
7 
8 static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC, Mat, Mat);
9 static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC, Mat, Mat);
10 static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC);
11 
12 static PetscErrorCode PCSetUp_BJacobi(PC pc) {
13   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
14   Mat         mat = pc->mat, pmat = pc->pmat;
15   PetscBool   hasop;
16   PetscInt    N, M, start, i, sum, end;
17   PetscInt    bs, i_start = -1, i_end = -1;
18   PetscMPIInt rank, size;
19 
20   PetscFunctionBegin;
21   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
22   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
23   PetscCall(MatGetLocalSize(pc->pmat, &M, &N));
24   PetscCall(MatGetBlockSize(pc->pmat, &bs));
25 
26   if (jac->n > 0 && jac->n < size) {
27     PetscCall(PCSetUp_BJacobi_Multiproc(pc));
28     PetscFunctionReturn(0);
29   }
30 
31   /*    Determines the number of blocks assigned to each processor */
32   /*   local block count  given */
33   if (jac->n_local > 0 && jac->n < 0) {
34     PetscCall(MPIU_Allreduce(&jac->n_local, &jac->n, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)pc)));
35     if (jac->l_lens) { /* check that user set these correctly */
36       sum = 0;
37       for (i = 0; i < jac->n_local; i++) {
38         PetscCheck(jac->l_lens[i] / bs * bs == jac->l_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat blocksize doesn't match block Jacobi layout");
39         sum += jac->l_lens[i];
40       }
41       PetscCheck(sum == M, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local lens set incorrectly");
42     } else {
43       PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
44       for (i = 0; i < jac->n_local; i++) jac->l_lens[i] = bs * ((M / bs) / jac->n_local + (((M / bs) % jac->n_local) > i));
45     }
46   } else if (jac->n > 0 && jac->n_local < 0) { /* global block count given */
47     /* global blocks given: determine which ones are local */
48     if (jac->g_lens) {
49       /* check if the g_lens is has valid entries */
50       for (i = 0; i < jac->n; i++) {
51         PetscCheck(jac->g_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Zero block not allowed");
52         PetscCheck(jac->g_lens[i] / bs * bs == jac->g_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat blocksize doesn't match block Jacobi layout");
53       }
54       if (size == 1) {
55         jac->n_local = jac->n;
56         PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
57         PetscCall(PetscArraycpy(jac->l_lens, jac->g_lens, jac->n_local));
58         /* check that user set these correctly */
59         sum = 0;
60         for (i = 0; i < jac->n_local; i++) sum += jac->l_lens[i];
61         PetscCheck(sum == M, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Global lens set incorrectly");
62       } else {
63         PetscCall(MatGetOwnershipRange(pc->pmat, &start, &end));
64         /* loop over blocks determing first one owned by me */
65         sum = 0;
66         for (i = 0; i < jac->n + 1; i++) {
67           if (sum == start) {
68             i_start = i;
69             goto start_1;
70           }
71           if (i < jac->n) sum += jac->g_lens[i];
72         }
73         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Block sizes used in PCBJacobiSetTotalBlocks()\nare not compatible with parallel matrix layout");
74       start_1:
75         for (i = i_start; i < jac->n + 1; i++) {
76           if (sum == end) {
77             i_end = i;
78             goto end_1;
79           }
80           if (i < jac->n) sum += jac->g_lens[i];
81         }
82         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Block sizes used in PCBJacobiSetTotalBlocks()\nare not compatible with parallel matrix layout");
83       end_1:
84         jac->n_local = i_end - i_start;
85         PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
86         PetscCall(PetscArraycpy(jac->l_lens, jac->g_lens + i_start, jac->n_local));
87       }
88     } else { /* no global blocks given, determine then using default layout */
89       jac->n_local = jac->n / size + ((jac->n % size) > rank);
90       PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
91       for (i = 0; i < jac->n_local; i++) {
92         jac->l_lens[i] = ((M / bs) / jac->n_local + (((M / bs) % jac->n_local) > i)) * bs;
93         PetscCheck(jac->l_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Too many blocks given");
94       }
95     }
96   } else if (jac->n < 0 && jac->n_local < 0) { /* no blocks given */
97     jac->n       = size;
98     jac->n_local = 1;
99     PetscCall(PetscMalloc1(1, &jac->l_lens));
100     jac->l_lens[0] = M;
101   } else { /* jac->n > 0 && jac->n_local > 0 */
102     if (!jac->l_lens) {
103       PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
104       for (i = 0; i < jac->n_local; i++) jac->l_lens[i] = bs * ((M / bs) / jac->n_local + (((M / bs) % jac->n_local) > i));
105     }
106   }
107   PetscCheck(jac->n_local >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of blocks is less than number of processors");
108 
109   /*    Determines mat and pmat */
110   PetscCall(MatHasOperation(pc->mat, MATOP_GET_DIAGONAL_BLOCK, &hasop));
111   if (!hasop && size == 1) {
112     mat  = pc->mat;
113     pmat = pc->pmat;
114   } else {
115     if (pc->useAmat) {
116       /* use block from Amat matrix, not Pmat for local MatMult() */
117       PetscCall(MatGetDiagonalBlock(pc->mat, &mat));
118     }
119     if (pc->pmat != pc->mat || !pc->useAmat) {
120       PetscCall(MatGetDiagonalBlock(pc->pmat, &pmat));
121     } else pmat = mat;
122   }
123 
124   /*
125      Setup code depends on the number of blocks
126   */
127   if (jac->n_local == 1) {
128     PetscCall(PCSetUp_BJacobi_Singleblock(pc, mat, pmat));
129   } else {
130     PetscCall(PCSetUp_BJacobi_Multiblock(pc, mat, pmat));
131   }
132   PetscFunctionReturn(0);
133 }
134 
135 /* Default destroy, if it has never been setup */
136 static PetscErrorCode PCDestroy_BJacobi(PC pc) {
137   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
138 
139   PetscFunctionBegin;
140   PetscCall(PetscFree(jac->g_lens));
141   PetscCall(PetscFree(jac->l_lens));
142   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetSubKSP_C", NULL));
143   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetTotalBlocks_C", NULL));
144   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetTotalBlocks_C", NULL));
145   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetLocalBlocks_C", NULL));
146   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetLocalBlocks_C", NULL));
147   PetscCall(PetscFree(pc->data));
148   PetscFunctionReturn(0);
149 }
150 
151 static PetscErrorCode PCSetFromOptions_BJacobi(PC pc, PetscOptionItems *PetscOptionsObject) {
152   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
153   PetscInt    blocks, i;
154   PetscBool   flg;
155 
156   PetscFunctionBegin;
157   PetscOptionsHeadBegin(PetscOptionsObject, "Block Jacobi options");
158   PetscCall(PetscOptionsInt("-pc_bjacobi_blocks", "Total number of blocks", "PCBJacobiSetTotalBlocks", jac->n, &blocks, &flg));
159   if (flg) PetscCall(PCBJacobiSetTotalBlocks(pc, blocks, NULL));
160   PetscCall(PetscOptionsInt("-pc_bjacobi_local_blocks", "Local number of blocks", "PCBJacobiSetLocalBlocks", jac->n_local, &blocks, &flg));
161   if (flg) PetscCall(PCBJacobiSetLocalBlocks(pc, blocks, NULL));
162   if (jac->ksp) {
163     /* The sub-KSP has already been set up (e.g., PCSetUp_BJacobi_Singleblock), but KSPSetFromOptions was not called
164      * unless we had already been called. */
165     for (i = 0; i < jac->n_local; i++) PetscCall(KSPSetFromOptions(jac->ksp[i]));
166   }
167   PetscOptionsHeadEnd();
168   PetscFunctionReturn(0);
169 }
170 
171 #include <petscdraw.h>
172 static PetscErrorCode PCView_BJacobi(PC pc, PetscViewer viewer) {
173   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
174   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
175   PetscMPIInt           rank;
176   PetscInt              i;
177   PetscBool             iascii, isstring, isdraw;
178   PetscViewer           sviewer;
179   PetscViewerFormat     format;
180   const char           *prefix;
181 
182   PetscFunctionBegin;
183   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
184   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
185   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
186   if (iascii) {
187     if (pc->useAmat) PetscCall(PetscViewerASCIIPrintf(viewer, "  using Amat local matrix, number of blocks = %" PetscInt_FMT "\n", jac->n));
188     PetscCall(PetscViewerASCIIPrintf(viewer, "  number of blocks = %" PetscInt_FMT "\n", jac->n));
189     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
190     PetscCall(PetscViewerGetFormat(viewer, &format));
191     if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) {
192       PetscCall(PetscViewerASCIIPrintf(viewer, "  Local solver information for first block is in the following KSP and PC objects on rank 0:\n"));
193       PetscCall(PCGetOptionsPrefix(pc, &prefix));
194       PetscCall(PetscViewerASCIIPrintf(viewer, "  Use -%sksp_view ::ascii_info_detail to display information for all blocks\n", prefix ? prefix : ""));
195       if (jac->ksp && !jac->psubcomm) {
196         PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
197         if (rank == 0) {
198           PetscCall(PetscViewerASCIIPushTab(viewer));
199           PetscCall(KSPView(jac->ksp[0], sviewer));
200           PetscCall(PetscViewerASCIIPopTab(viewer));
201         }
202         PetscCall(PetscViewerFlush(sviewer));
203         PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
204         PetscCall(PetscViewerFlush(viewer));
205         /*  extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
206         PetscCall(PetscViewerASCIIPopSynchronized(viewer));
207       } else if (mpjac && jac->ksp && mpjac->psubcomm) {
208         PetscCall(PetscViewerGetSubViewer(viewer, mpjac->psubcomm->child, &sviewer));
209         if (!mpjac->psubcomm->color) {
210           PetscCall(PetscViewerASCIIPushTab(viewer));
211           PetscCall(KSPView(*(jac->ksp), sviewer));
212           PetscCall(PetscViewerASCIIPopTab(viewer));
213         }
214         PetscCall(PetscViewerFlush(sviewer));
215         PetscCall(PetscViewerRestoreSubViewer(viewer, mpjac->psubcomm->child, &sviewer));
216         PetscCall(PetscViewerFlush(viewer));
217         /*  extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
218         PetscCall(PetscViewerASCIIPopSynchronized(viewer));
219       } else {
220         PetscCall(PetscViewerFlush(viewer));
221       }
222     } else {
223       PetscInt n_global;
224       PetscCall(MPIU_Allreduce(&jac->n_local, &n_global, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)pc)));
225       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
226       PetscCall(PetscViewerASCIIPrintf(viewer, "  Local solver information for each block is in the following KSP and PC objects:\n"));
227       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] number of local blocks = %" PetscInt_FMT ", first local block number = %" PetscInt_FMT "\n", rank, jac->n_local, jac->first_local));
228       PetscCall(PetscViewerASCIIPushTab(viewer));
229       PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
230       for (i = 0; i < jac->n_local; i++) {
231         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] local block number %" PetscInt_FMT "\n", rank, i));
232         PetscCall(KSPView(jac->ksp[i], sviewer));
233         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "- - - - - - - - - - - - - - - - - -\n"));
234       }
235       PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
236       PetscCall(PetscViewerASCIIPopTab(viewer));
237       PetscCall(PetscViewerFlush(viewer));
238       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
239     }
240   } else if (isstring) {
241     PetscCall(PetscViewerStringSPrintf(viewer, " blks=%" PetscInt_FMT, jac->n));
242     PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
243     if (jac->ksp) PetscCall(KSPView(jac->ksp[0], sviewer));
244     PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
245   } else if (isdraw) {
246     PetscDraw draw;
247     char      str[25];
248     PetscReal x, y, bottom, h;
249 
250     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
251     PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
252     PetscCall(PetscSNPrintf(str, 25, "Number blocks %" PetscInt_FMT, jac->n));
253     PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h));
254     bottom = y - h;
255     PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
256     /* warning the communicator on viewer is different then on ksp in parallel */
257     if (jac->ksp) PetscCall(KSPView(jac->ksp[0], viewer));
258     PetscCall(PetscDrawPopCurrentPoint(draw));
259   }
260   PetscFunctionReturn(0);
261 }
262 
263 static PetscErrorCode PCBJacobiGetSubKSP_BJacobi(PC pc, PetscInt *n_local, PetscInt *first_local, KSP **ksp) {
264   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
265 
266   PetscFunctionBegin;
267   PetscCheck(pc->setupcalled, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must call KSPSetUp() or PCSetUp() first");
268 
269   if (n_local) *n_local = jac->n_local;
270   if (first_local) *first_local = jac->first_local;
271   if (ksp) *ksp = jac->ksp;
272   PetscFunctionReturn(0);
273 }
274 
275 static PetscErrorCode PCBJacobiSetTotalBlocks_BJacobi(PC pc, PetscInt blocks, PetscInt *lens) {
276   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
277 
278   PetscFunctionBegin;
279   PetscCheck(pc->setupcalled <= 0 || jac->n == blocks, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Cannot alter number of blocks after PCSetUp()/KSPSetUp() has been called");
280   jac->n = blocks;
281   if (!lens) jac->g_lens = NULL;
282   else {
283     PetscCall(PetscMalloc1(blocks, &jac->g_lens));
284     PetscCall(PetscLogObjectMemory((PetscObject)pc, blocks * sizeof(PetscInt)));
285     PetscCall(PetscArraycpy(jac->g_lens, lens, blocks));
286   }
287   PetscFunctionReturn(0);
288 }
289 
290 static PetscErrorCode PCBJacobiGetTotalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[]) {
291   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
292 
293   PetscFunctionBegin;
294   *blocks = jac->n;
295   if (lens) *lens = jac->g_lens;
296   PetscFunctionReturn(0);
297 }
298 
299 static PetscErrorCode PCBJacobiSetLocalBlocks_BJacobi(PC pc, PetscInt blocks, const PetscInt lens[]) {
300   PC_BJacobi *jac;
301 
302   PetscFunctionBegin;
303   jac = (PC_BJacobi *)pc->data;
304 
305   jac->n_local = blocks;
306   if (!lens) jac->l_lens = NULL;
307   else {
308     PetscCall(PetscMalloc1(blocks, &jac->l_lens));
309     PetscCall(PetscLogObjectMemory((PetscObject)pc, blocks * sizeof(PetscInt)));
310     PetscCall(PetscArraycpy(jac->l_lens, lens, blocks));
311   }
312   PetscFunctionReturn(0);
313 }
314 
315 static PetscErrorCode PCBJacobiGetLocalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[]) {
316   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
317 
318   PetscFunctionBegin;
319   *blocks = jac->n_local;
320   if (lens) *lens = jac->l_lens;
321   PetscFunctionReturn(0);
322 }
323 
324 /*@C
325    PCBJacobiGetSubKSP - Gets the local `KSP` contexts for all blocks on
326    this processor.
327 
328    Not Collective
329 
330    Input Parameter:
331 .  pc - the preconditioner context
332 
333    Output Parameters:
334 +  n_local - the number of blocks on this processor, or NULL
335 .  first_local - the global number of the first block on this processor, or NULL
336 -  ksp - the array of KSP contexts
337 
338    Notes:
339    After `PCBJacobiGetSubKSP()` the array of `KSP` contexts is not to be freed.
340 
341    Currently for some matrix implementations only 1 block per processor
342    is supported.
343 
344    You must call `KSPSetUp()` or `PCSetUp()` before calling `PCBJacobiGetSubKSP()`.
345 
346    Fortran Usage:
347    You must pass in a `KSP` array that is large enough to contain all the local `KSP`s.
348 
349    You can call `PCBJacobiGetSubKSP`(pc,nlocal,firstlocal,`PETSC_NULL_KSP`,ierr) to determine how large the
350    `KSP` array must be.
351 
352    Level: advanced
353 
354 .seealso: `PCBJACOBI`, `PCASM`, `PCASMGetSubKSP()`
355 @*/
356 PetscErrorCode PCBJacobiGetSubKSP(PC pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp[]) {
357   PetscFunctionBegin;
358   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
359   PetscUseMethod(pc, "PCBJacobiGetSubKSP_C", (PC, PetscInt *, PetscInt *, KSP **), (pc, n_local, first_local, ksp));
360   PetscFunctionReturn(0);
361 }
362 
363 /*@
364    PCBJacobiSetTotalBlocks - Sets the global number of blocks for the block
365    Jacobi preconditioner.
366 
367    Collective on pc
368 
369    Input Parameters:
370 +  pc - the preconditioner context
371 .  blocks - the number of blocks
372 -  lens - [optional] integer array containing the size of each block
373 
374    Options Database Key:
375 .  -pc_bjacobi_blocks <blocks> - Sets the number of global blocks
376 
377    Note:
378    Currently only a limited number of blocking configurations are supported.
379    All processors sharing the `PC` must call this routine with the same data.
380 
381    Level: intermediate
382 
383 .seealso: `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiSetLocalBlocks()`
384 @*/
385 PetscErrorCode PCBJacobiSetTotalBlocks(PC pc, PetscInt blocks, const PetscInt lens[]) {
386   PetscFunctionBegin;
387   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
388   PetscCheck(blocks > 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Must have positive blocks");
389   PetscTryMethod(pc, "PCBJacobiSetTotalBlocks_C", (PC, PetscInt, const PetscInt[]), (pc, blocks, lens));
390   PetscFunctionReturn(0);
391 }
392 
393 /*@C
394    PCBJacobiGetTotalBlocks - Gets the global number of blocks for the block
395    Jacobi, `PCBJACOBI`, preconditioner.
396 
397    Not Collective
398 
399    Input Parameter:
400 .  pc - the preconditioner context
401 
402    Output parameters:
403 +  blocks - the number of blocks
404 -  lens - integer array containing the size of each block
405 
406    Level: intermediate
407 
408 .seealso: `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiGetLocalBlocks()`
409 @*/
410 PetscErrorCode PCBJacobiGetTotalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[]) {
411   PetscFunctionBegin;
412   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
413   PetscValidIntPointer(blocks, 2);
414   PetscUseMethod(pc, "PCBJacobiGetTotalBlocks_C", (PC, PetscInt *, const PetscInt *[]), (pc, blocks, lens));
415   PetscFunctionReturn(0);
416 }
417 
418 /*@
419    PCBJacobiSetLocalBlocks - Sets the local number of blocks for the block
420    Jacobi, `PCBJACOBI`,  preconditioner.
421 
422    Not Collective
423 
424    Input Parameters:
425 +  pc - the preconditioner context
426 .  blocks - the number of blocks
427 -  lens - [optional] integer array containing size of each block
428 
429    Options Database Key:
430 .  -pc_bjacobi_local_blocks <blocks> - Sets the number of local blocks
431 
432    Note:
433    Currently only a limited number of blocking configurations are supported.
434 
435    Level: intermediate
436 
437 .seealso: `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiSetTotalBlocks()`
438 @*/
439 PetscErrorCode PCBJacobiSetLocalBlocks(PC pc, PetscInt blocks, const PetscInt lens[]) {
440   PetscFunctionBegin;
441   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
442   PetscCheck(blocks >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have nonegative blocks");
443   PetscTryMethod(pc, "PCBJacobiSetLocalBlocks_C", (PC, PetscInt, const PetscInt[]), (pc, blocks, lens));
444   PetscFunctionReturn(0);
445 }
446 
447 /*@C
448    PCBJacobiGetLocalBlocks - Gets the local number of blocks for the block
449    Jacobi, `PCBJACOBI`, preconditioner.
450 
451    Not Collective
452 
453    Input Parameters:
454 +  pc - the preconditioner context
455 .  blocks - the number of blocks
456 -  lens - [optional] integer array containing size of each block
457 
458    Note:
459    Currently only a limited number of blocking configurations are supported.
460 
461    Level: intermediate
462 
463 .seealso: `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiGetTotalBlocks()`
464 @*/
465 PetscErrorCode PCBJacobiGetLocalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[]) {
466   PetscFunctionBegin;
467   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
468   PetscValidIntPointer(blocks, 2);
469   PetscUseMethod(pc, "PCBJacobiGetLocalBlocks_C", (PC, PetscInt *, const PetscInt *[]), (pc, blocks, lens));
470   PetscFunctionReturn(0);
471 }
472 
473 /*MC
474    PCBJACOBI - Use block Jacobi preconditioning, each block is (approximately) solved with
475            its own `KSP` object.
476 
477    Options Database Keys:
478 +  -pc_use_amat - use Amat to apply block of operator in inner Krylov method
479 -  -pc_bjacobi_blocks <n> - use n total blocks
480 
481    Notes:
482     See `PCJACOBI` for diagonal Jacobi, `PCVPBJACOBI` for variable point block, and `PCPBJACOBI` for fixed size point block
483 
484     Each processor can have one or more blocks, or a single block can be shared by several processes. Defaults to one block per processor.
485 
486      To set options on the solvers for each block append -sub_ to all the `KSP` and `PC`
487         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
488 
489      To set the options on the solvers separate for each block call `PCBJacobiGetSubKSP()`
490          and set the options directly on the resulting `KSP` object (you can access its `PC`
491          `KSPGetPC())`
492 
493      For GPU-based vectors (`VECCUDA`, `VECViennaCL`) it is recommended to use exactly one block per MPI process for best
494          performance.  Different block partitioning may lead to additional data transfers
495          between host and GPU that lead to degraded performance.
496 
497      When multiple processes share a single block, each block encompasses exactly all the unknowns owned its set of processes.
498 
499    Level: beginner
500 
501 .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCType`,
502           `PCASM`, `PCSetUseAmat()`, `PCGetUseAmat()`, `PCBJacobiGetSubKSP()`, `PCBJacobiSetTotalBlocks()`,
503           `PCBJacobiSetLocalBlocks()`, `PCSetModifySubMatrices()`, `PCJACOBI`, `PCVPBJACOBI`, `PCPBJACOBI`
504 M*/
505 
506 PETSC_EXTERN PetscErrorCode PCCreate_BJacobi(PC pc) {
507   PetscMPIInt rank;
508   PC_BJacobi *jac;
509 
510   PetscFunctionBegin;
511   PetscCall(PetscNewLog(pc, &jac));
512   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
513 
514   pc->ops->apply           = NULL;
515   pc->ops->matapply        = NULL;
516   pc->ops->applytranspose  = NULL;
517   pc->ops->setup           = PCSetUp_BJacobi;
518   pc->ops->destroy         = PCDestroy_BJacobi;
519   pc->ops->setfromoptions  = PCSetFromOptions_BJacobi;
520   pc->ops->view            = PCView_BJacobi;
521   pc->ops->applyrichardson = NULL;
522 
523   pc->data         = (void *)jac;
524   jac->n           = -1;
525   jac->n_local     = -1;
526   jac->first_local = rank;
527   jac->ksp         = NULL;
528   jac->g_lens      = NULL;
529   jac->l_lens      = NULL;
530   jac->psubcomm    = NULL;
531 
532   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetSubKSP_C", PCBJacobiGetSubKSP_BJacobi));
533   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetTotalBlocks_C", PCBJacobiSetTotalBlocks_BJacobi));
534   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetTotalBlocks_C", PCBJacobiGetTotalBlocks_BJacobi));
535   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetLocalBlocks_C", PCBJacobiSetLocalBlocks_BJacobi));
536   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetLocalBlocks_C", PCBJacobiGetLocalBlocks_BJacobi));
537   PetscFunctionReturn(0);
538 }
539 
540 /*
541         These are for a single block per processor; works for AIJ, BAIJ; Seq and MPI
542 */
543 static PetscErrorCode PCReset_BJacobi_Singleblock(PC pc) {
544   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
545   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
546 
547   PetscFunctionBegin;
548   PetscCall(KSPReset(jac->ksp[0]));
549   PetscCall(VecDestroy(&bjac->x));
550   PetscCall(VecDestroy(&bjac->y));
551   PetscFunctionReturn(0);
552 }
553 
554 static PetscErrorCode PCDestroy_BJacobi_Singleblock(PC pc) {
555   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
556   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
557 
558   PetscFunctionBegin;
559   PetscCall(PCReset_BJacobi_Singleblock(pc));
560   PetscCall(KSPDestroy(&jac->ksp[0]));
561   PetscCall(PetscFree(jac->ksp));
562   PetscCall(PetscFree(bjac));
563   PetscCall(PCDestroy_BJacobi(pc));
564   PetscFunctionReturn(0);
565 }
566 
567 static PetscErrorCode PCSetUpOnBlocks_BJacobi_Singleblock(PC pc) {
568   PC_BJacobi        *jac    = (PC_BJacobi *)pc->data;
569   KSP                subksp = jac->ksp[0];
570   KSPConvergedReason reason;
571 
572   PetscFunctionBegin;
573   PetscCall(KSPSetUp(subksp));
574   PetscCall(KSPGetConvergedReason(subksp, &reason));
575   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
576   PetscFunctionReturn(0);
577 }
578 
579 static PetscErrorCode PCApply_BJacobi_Singleblock(PC pc, Vec x, Vec y) {
580   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
581   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
582 
583   PetscFunctionBegin;
584   PetscCall(VecGetLocalVectorRead(x, bjac->x));
585   PetscCall(VecGetLocalVector(y, bjac->y));
586   /* Since the inner KSP matrix may point directly to the diagonal block of an MPI matrix the inner
587      matrix may change even if the outer KSP/PC has not updated the preconditioner, this will trigger a rebuild
588      of the inner preconditioner automatically unless we pass down the outer preconditioners reuse flag.*/
589   PetscCall(KSPSetReusePreconditioner(jac->ksp[0], pc->reusepreconditioner));
590   PetscCall(KSPSolve(jac->ksp[0], bjac->x, bjac->y));
591   PetscCall(KSPCheckSolve(jac->ksp[0], pc, bjac->y));
592   PetscCall(VecRestoreLocalVectorRead(x, bjac->x));
593   PetscCall(VecRestoreLocalVector(y, bjac->y));
594   PetscFunctionReturn(0);
595 }
596 
597 static PetscErrorCode PCMatApply_BJacobi_Singleblock(PC pc, Mat X, Mat Y) {
598   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
599   Mat         sX, sY;
600 
601   PetscFunctionBegin;
602   /* Since the inner KSP matrix may point directly to the diagonal block of an MPI matrix the inner
603      matrix may change even if the outer KSP/PC has not updated the preconditioner, this will trigger a rebuild
604      of the inner preconditioner automatically unless we pass down the outer preconditioners reuse flag.*/
605   PetscCall(KSPSetReusePreconditioner(jac->ksp[0], pc->reusepreconditioner));
606   PetscCall(MatDenseGetLocalMatrix(X, &sX));
607   PetscCall(MatDenseGetLocalMatrix(Y, &sY));
608   PetscCall(KSPMatSolve(jac->ksp[0], sX, sY));
609   PetscFunctionReturn(0);
610 }
611 
612 static PetscErrorCode PCApplySymmetricLeft_BJacobi_Singleblock(PC pc, Vec x, Vec y) {
613   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
614   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
615   PetscScalar            *y_array;
616   const PetscScalar      *x_array;
617   PC                      subpc;
618 
619   PetscFunctionBegin;
620   /*
621       The VecPlaceArray() is to avoid having to copy the
622     y vector into the bjac->x vector. The reason for
623     the bjac->x vector is that we need a sequential vector
624     for the sequential solve.
625   */
626   PetscCall(VecGetArrayRead(x, &x_array));
627   PetscCall(VecGetArray(y, &y_array));
628   PetscCall(VecPlaceArray(bjac->x, x_array));
629   PetscCall(VecPlaceArray(bjac->y, y_array));
630   /* apply the symmetric left portion of the inner PC operator */
631   /* note this by-passes the inner KSP and its options completely */
632   PetscCall(KSPGetPC(jac->ksp[0], &subpc));
633   PetscCall(PCApplySymmetricLeft(subpc, bjac->x, bjac->y));
634   PetscCall(VecResetArray(bjac->x));
635   PetscCall(VecResetArray(bjac->y));
636   PetscCall(VecRestoreArrayRead(x, &x_array));
637   PetscCall(VecRestoreArray(y, &y_array));
638   PetscFunctionReturn(0);
639 }
640 
641 static PetscErrorCode PCApplySymmetricRight_BJacobi_Singleblock(PC pc, Vec x, Vec y) {
642   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
643   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
644   PetscScalar            *y_array;
645   const PetscScalar      *x_array;
646   PC                      subpc;
647 
648   PetscFunctionBegin;
649   /*
650       The VecPlaceArray() is to avoid having to copy the
651     y vector into the bjac->x vector. The reason for
652     the bjac->x vector is that we need a sequential vector
653     for the sequential solve.
654   */
655   PetscCall(VecGetArrayRead(x, &x_array));
656   PetscCall(VecGetArray(y, &y_array));
657   PetscCall(VecPlaceArray(bjac->x, x_array));
658   PetscCall(VecPlaceArray(bjac->y, y_array));
659 
660   /* apply the symmetric right portion of the inner PC operator */
661   /* note this by-passes the inner KSP and its options completely */
662 
663   PetscCall(KSPGetPC(jac->ksp[0], &subpc));
664   PetscCall(PCApplySymmetricRight(subpc, bjac->x, bjac->y));
665 
666   PetscCall(VecRestoreArrayRead(x, &x_array));
667   PetscCall(VecRestoreArray(y, &y_array));
668   PetscFunctionReturn(0);
669 }
670 
671 static PetscErrorCode PCApplyTranspose_BJacobi_Singleblock(PC pc, Vec x, Vec y) {
672   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
673   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
674   PetscScalar            *y_array;
675   const PetscScalar      *x_array;
676 
677   PetscFunctionBegin;
678   /*
679       The VecPlaceArray() is to avoid having to copy the
680     y vector into the bjac->x vector. The reason for
681     the bjac->x vector is that we need a sequential vector
682     for the sequential solve.
683   */
684   PetscCall(VecGetArrayRead(x, &x_array));
685   PetscCall(VecGetArray(y, &y_array));
686   PetscCall(VecPlaceArray(bjac->x, x_array));
687   PetscCall(VecPlaceArray(bjac->y, y_array));
688   PetscCall(KSPSolveTranspose(jac->ksp[0], bjac->x, bjac->y));
689   PetscCall(KSPCheckSolve(jac->ksp[0], pc, bjac->y));
690   PetscCall(VecResetArray(bjac->x));
691   PetscCall(VecResetArray(bjac->y));
692   PetscCall(VecRestoreArrayRead(x, &x_array));
693   PetscCall(VecRestoreArray(y, &y_array));
694   PetscFunctionReturn(0);
695 }
696 
697 static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC pc, Mat mat, Mat pmat) {
698   PC_BJacobi             *jac = (PC_BJacobi *)pc->data;
699   PetscInt                m;
700   KSP                     ksp;
701   PC_BJacobi_Singleblock *bjac;
702   PetscBool               wasSetup = PETSC_TRUE;
703   VecType                 vectype;
704   const char             *prefix;
705 
706   PetscFunctionBegin;
707   if (!pc->setupcalled) {
708     if (!jac->ksp) {
709       wasSetup = PETSC_FALSE;
710 
711       PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
712       PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure));
713       PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1));
714       PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)ksp));
715       PetscCall(KSPSetType(ksp, KSPPREONLY));
716       PetscCall(PCGetOptionsPrefix(pc, &prefix));
717       PetscCall(KSPSetOptionsPrefix(ksp, prefix));
718       PetscCall(KSPAppendOptionsPrefix(ksp, "sub_"));
719 
720       pc->ops->reset               = PCReset_BJacobi_Singleblock;
721       pc->ops->destroy             = PCDestroy_BJacobi_Singleblock;
722       pc->ops->apply               = PCApply_BJacobi_Singleblock;
723       pc->ops->matapply            = PCMatApply_BJacobi_Singleblock;
724       pc->ops->applysymmetricleft  = PCApplySymmetricLeft_BJacobi_Singleblock;
725       pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Singleblock;
726       pc->ops->applytranspose      = PCApplyTranspose_BJacobi_Singleblock;
727       pc->ops->setuponblocks       = PCSetUpOnBlocks_BJacobi_Singleblock;
728 
729       PetscCall(PetscMalloc1(1, &jac->ksp));
730       jac->ksp[0] = ksp;
731 
732       PetscCall(PetscNewLog(pc, &bjac));
733       jac->data = (void *)bjac;
734     } else {
735       ksp  = jac->ksp[0];
736       bjac = (PC_BJacobi_Singleblock *)jac->data;
737     }
738 
739     /*
740       The reason we need to generate these vectors is to serve
741       as the right-hand side and solution vector for the solve on the
742       block. We do not need to allocate space for the vectors since
743       that is provided via VecPlaceArray() just before the call to
744       KSPSolve() on the block.
745     */
746     PetscCall(MatGetSize(pmat, &m, &m));
747     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &bjac->x));
748     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &bjac->y));
749     PetscCall(MatGetVecType(pmat, &vectype));
750     PetscCall(VecSetType(bjac->x, vectype));
751     PetscCall(VecSetType(bjac->y, vectype));
752     PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)bjac->x));
753     PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)bjac->y));
754   } else {
755     ksp  = jac->ksp[0];
756     bjac = (PC_BJacobi_Singleblock *)jac->data;
757   }
758   PetscCall(KSPGetOptionsPrefix(ksp, &prefix));
759   if (pc->useAmat) {
760     PetscCall(KSPSetOperators(ksp, mat, pmat));
761     PetscCall(MatSetOptionsPrefix(mat, prefix));
762   } else {
763     PetscCall(KSPSetOperators(ksp, pmat, pmat));
764   }
765   PetscCall(MatSetOptionsPrefix(pmat, prefix));
766   if (!wasSetup && pc->setfromoptionscalled) {
767     /* If PCSetFromOptions_BJacobi is called later, KSPSetFromOptions will be called at that time. */
768     PetscCall(KSPSetFromOptions(ksp));
769   }
770   PetscFunctionReturn(0);
771 }
772 
773 static PetscErrorCode PCReset_BJacobi_Multiblock(PC pc) {
774   PC_BJacobi            *jac  = (PC_BJacobi *)pc->data;
775   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
776   PetscInt               i;
777 
778   PetscFunctionBegin;
779   if (bjac && bjac->pmat) {
780     PetscCall(MatDestroyMatrices(jac->n_local, &bjac->pmat));
781     if (pc->useAmat) PetscCall(MatDestroyMatrices(jac->n_local, &bjac->mat));
782   }
783 
784   for (i = 0; i < jac->n_local; i++) {
785     PetscCall(KSPReset(jac->ksp[i]));
786     if (bjac && bjac->x) {
787       PetscCall(VecDestroy(&bjac->x[i]));
788       PetscCall(VecDestroy(&bjac->y[i]));
789       PetscCall(ISDestroy(&bjac->is[i]));
790     }
791   }
792   PetscCall(PetscFree(jac->l_lens));
793   PetscCall(PetscFree(jac->g_lens));
794   PetscFunctionReturn(0);
795 }
796 
797 static PetscErrorCode PCDestroy_BJacobi_Multiblock(PC pc) {
798   PC_BJacobi            *jac  = (PC_BJacobi *)pc->data;
799   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
800   PetscInt               i;
801 
802   PetscFunctionBegin;
803   PetscCall(PCReset_BJacobi_Multiblock(pc));
804   if (bjac) {
805     PetscCall(PetscFree2(bjac->x, bjac->y));
806     PetscCall(PetscFree(bjac->starts));
807     PetscCall(PetscFree(bjac->is));
808   }
809   PetscCall(PetscFree(jac->data));
810   for (i = 0; i < jac->n_local; i++) PetscCall(KSPDestroy(&jac->ksp[i]));
811   PetscCall(PetscFree(jac->ksp));
812   PetscCall(PCDestroy_BJacobi(pc));
813   PetscFunctionReturn(0);
814 }
815 
816 static PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiblock(PC pc) {
817   PC_BJacobi        *jac = (PC_BJacobi *)pc->data;
818   PetscInt           i, n_local = jac->n_local;
819   KSPConvergedReason reason;
820 
821   PetscFunctionBegin;
822   for (i = 0; i < n_local; i++) {
823     PetscCall(KSPSetUp(jac->ksp[i]));
824     PetscCall(KSPGetConvergedReason(jac->ksp[i], &reason));
825     if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
826   }
827   PetscFunctionReturn(0);
828 }
829 
830 /*
831       Preconditioner for block Jacobi
832 */
833 static PetscErrorCode PCApply_BJacobi_Multiblock(PC pc, Vec x, Vec y) {
834   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
835   PetscInt               i, n_local = jac->n_local;
836   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
837   PetscScalar           *yin;
838   const PetscScalar     *xin;
839 
840   PetscFunctionBegin;
841   PetscCall(VecGetArrayRead(x, &xin));
842   PetscCall(VecGetArray(y, &yin));
843   for (i = 0; i < n_local; i++) {
844     /*
845        To avoid copying the subvector from x into a workspace we instead
846        make the workspace vector array point to the subpart of the array of
847        the global vector.
848     */
849     PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
850     PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
851 
852     PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
853     PetscCall(KSPSolve(jac->ksp[i], bjac->x[i], bjac->y[i]));
854     PetscCall(KSPCheckSolve(jac->ksp[i], pc, bjac->y[i]));
855     PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
856 
857     PetscCall(VecResetArray(bjac->x[i]));
858     PetscCall(VecResetArray(bjac->y[i]));
859   }
860   PetscCall(VecRestoreArrayRead(x, &xin));
861   PetscCall(VecRestoreArray(y, &yin));
862   PetscFunctionReturn(0);
863 }
864 
865 /*
866       Preconditioner for block Jacobi
867 */
868 static PetscErrorCode PCApplyTranspose_BJacobi_Multiblock(PC pc, Vec x, Vec y) {
869   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
870   PetscInt               i, n_local = jac->n_local;
871   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
872   PetscScalar           *yin;
873   const PetscScalar     *xin;
874 
875   PetscFunctionBegin;
876   PetscCall(VecGetArrayRead(x, &xin));
877   PetscCall(VecGetArray(y, &yin));
878   for (i = 0; i < n_local; i++) {
879     /*
880        To avoid copying the subvector from x into a workspace we instead
881        make the workspace vector array point to the subpart of the array of
882        the global vector.
883     */
884     PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
885     PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
886 
887     PetscCall(PetscLogEventBegin(PC_ApplyTransposeOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
888     PetscCall(KSPSolveTranspose(jac->ksp[i], bjac->x[i], bjac->y[i]));
889     PetscCall(KSPCheckSolve(jac->ksp[i], pc, bjac->y[i]));
890     PetscCall(PetscLogEventEnd(PC_ApplyTransposeOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
891 
892     PetscCall(VecResetArray(bjac->x[i]));
893     PetscCall(VecResetArray(bjac->y[i]));
894   }
895   PetscCall(VecRestoreArrayRead(x, &xin));
896   PetscCall(VecRestoreArray(y, &yin));
897   PetscFunctionReturn(0);
898 }
899 
900 static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC pc, Mat mat, Mat pmat) {
901   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
902   PetscInt               m, n_local, N, M, start, i;
903   const char            *prefix;
904   KSP                    ksp;
905   Vec                    x, y;
906   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
907   PC                     subpc;
908   IS                     is;
909   MatReuse               scall;
910   VecType                vectype;
911 
912   PetscFunctionBegin;
913   PetscCall(MatGetLocalSize(pc->pmat, &M, &N));
914 
915   n_local = jac->n_local;
916 
917   if (pc->useAmat) {
918     PetscBool same;
919     PetscCall(PetscObjectTypeCompare((PetscObject)mat, ((PetscObject)pmat)->type_name, &same));
920     PetscCheck(same, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Matrices not of same type");
921   }
922 
923   if (!pc->setupcalled) {
924     scall = MAT_INITIAL_MATRIX;
925 
926     if (!jac->ksp) {
927       pc->ops->reset          = PCReset_BJacobi_Multiblock;
928       pc->ops->destroy        = PCDestroy_BJacobi_Multiblock;
929       pc->ops->apply          = PCApply_BJacobi_Multiblock;
930       pc->ops->matapply       = NULL;
931       pc->ops->applytranspose = PCApplyTranspose_BJacobi_Multiblock;
932       pc->ops->setuponblocks  = PCSetUpOnBlocks_BJacobi_Multiblock;
933 
934       PetscCall(PetscNewLog(pc, &bjac));
935       PetscCall(PetscMalloc1(n_local, &jac->ksp));
936       PetscCall(PetscLogObjectMemory((PetscObject)pc, sizeof(n_local * sizeof(KSP))));
937       PetscCall(PetscMalloc2(n_local, &bjac->x, n_local, &bjac->y));
938       PetscCall(PetscMalloc1(n_local, &bjac->starts));
939       PetscCall(PetscLogObjectMemory((PetscObject)pc, sizeof(n_local * sizeof(PetscScalar))));
940 
941       jac->data = (void *)bjac;
942       PetscCall(PetscMalloc1(n_local, &bjac->is));
943       PetscCall(PetscLogObjectMemory((PetscObject)pc, sizeof(n_local * sizeof(IS))));
944 
945       for (i = 0; i < n_local; i++) {
946         PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
947         PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure));
948         PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1));
949         PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)ksp));
950         PetscCall(KSPSetType(ksp, KSPPREONLY));
951         PetscCall(KSPGetPC(ksp, &subpc));
952         PetscCall(PCGetOptionsPrefix(pc, &prefix));
953         PetscCall(KSPSetOptionsPrefix(ksp, prefix));
954         PetscCall(KSPAppendOptionsPrefix(ksp, "sub_"));
955 
956         jac->ksp[i] = ksp;
957       }
958     } else {
959       bjac = (PC_BJacobi_Multiblock *)jac->data;
960     }
961 
962     start = 0;
963     PetscCall(MatGetVecType(pmat, &vectype));
964     for (i = 0; i < n_local; i++) {
965       m = jac->l_lens[i];
966       /*
967       The reason we need to generate these vectors is to serve
968       as the right-hand side and solution vector for the solve on the
969       block. We do not need to allocate space for the vectors since
970       that is provided via VecPlaceArray() just before the call to
971       KSPSolve() on the block.
972 
973       */
974       PetscCall(VecCreateSeq(PETSC_COMM_SELF, m, &x));
975       PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &y));
976       PetscCall(VecSetType(x, vectype));
977       PetscCall(VecSetType(y, vectype));
978       PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)x));
979       PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)y));
980 
981       bjac->x[i]      = x;
982       bjac->y[i]      = y;
983       bjac->starts[i] = start;
984 
985       PetscCall(ISCreateStride(PETSC_COMM_SELF, m, start, 1, &is));
986       bjac->is[i] = is;
987       PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)is));
988 
989       start += m;
990     }
991   } else {
992     bjac = (PC_BJacobi_Multiblock *)jac->data;
993     /*
994        Destroy the blocks from the previous iteration
995     */
996     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
997       PetscCall(MatDestroyMatrices(n_local, &bjac->pmat));
998       if (pc->useAmat) PetscCall(MatDestroyMatrices(n_local, &bjac->mat));
999       scall = MAT_INITIAL_MATRIX;
1000     } else scall = MAT_REUSE_MATRIX;
1001   }
1002 
1003   PetscCall(MatCreateSubMatrices(pmat, n_local, bjac->is, bjac->is, scall, &bjac->pmat));
1004   if (pc->useAmat) PetscCall(MatCreateSubMatrices(mat, n_local, bjac->is, bjac->is, scall, &bjac->mat));
1005   /* Return control to the user so that the submatrices can be modified (e.g., to apply
1006      different boundary conditions for the submatrices than for the global problem) */
1007   PetscCall(PCModifySubMatrices(pc, n_local, bjac->is, bjac->is, bjac->pmat, pc->modifysubmatricesP));
1008 
1009   for (i = 0; i < n_local; i++) {
1010     PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)bjac->pmat[i]));
1011     PetscCall(KSPGetOptionsPrefix(jac->ksp[i], &prefix));
1012     if (pc->useAmat) {
1013       PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)bjac->mat[i]));
1014       PetscCall(KSPSetOperators(jac->ksp[i], bjac->mat[i], bjac->pmat[i]));
1015       PetscCall(MatSetOptionsPrefix(bjac->mat[i], prefix));
1016     } else {
1017       PetscCall(KSPSetOperators(jac->ksp[i], bjac->pmat[i], bjac->pmat[i]));
1018     }
1019     PetscCall(MatSetOptionsPrefix(bjac->pmat[i], prefix));
1020     if (pc->setfromoptionscalled) PetscCall(KSPSetFromOptions(jac->ksp[i]));
1021   }
1022   PetscFunctionReturn(0);
1023 }
1024 
1025 /*
1026       These are for a single block with multiple processes
1027 */
1028 static PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiproc(PC pc) {
1029   PC_BJacobi        *jac    = (PC_BJacobi *)pc->data;
1030   KSP                subksp = jac->ksp[0];
1031   KSPConvergedReason reason;
1032 
1033   PetscFunctionBegin;
1034   PetscCall(KSPSetUp(subksp));
1035   PetscCall(KSPGetConvergedReason(subksp, &reason));
1036   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
1037   PetscFunctionReturn(0);
1038 }
1039 
1040 static PetscErrorCode PCReset_BJacobi_Multiproc(PC pc) {
1041   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
1042   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
1043 
1044   PetscFunctionBegin;
1045   PetscCall(VecDestroy(&mpjac->ysub));
1046   PetscCall(VecDestroy(&mpjac->xsub));
1047   PetscCall(MatDestroy(&mpjac->submats));
1048   if (jac->ksp) PetscCall(KSPReset(jac->ksp[0]));
1049   PetscFunctionReturn(0);
1050 }
1051 
1052 static PetscErrorCode PCDestroy_BJacobi_Multiproc(PC pc) {
1053   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
1054   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
1055 
1056   PetscFunctionBegin;
1057   PetscCall(PCReset_BJacobi_Multiproc(pc));
1058   PetscCall(KSPDestroy(&jac->ksp[0]));
1059   PetscCall(PetscFree(jac->ksp));
1060   PetscCall(PetscSubcommDestroy(&mpjac->psubcomm));
1061 
1062   PetscCall(PetscFree(mpjac));
1063   PetscCall(PCDestroy_BJacobi(pc));
1064   PetscFunctionReturn(0);
1065 }
1066 
1067 static PetscErrorCode PCApply_BJacobi_Multiproc(PC pc, Vec x, Vec y) {
1068   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
1069   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
1070   PetscScalar          *yarray;
1071   const PetscScalar    *xarray;
1072   KSPConvergedReason    reason;
1073 
1074   PetscFunctionBegin;
1075   /* place x's and y's local arrays into xsub and ysub */
1076   PetscCall(VecGetArrayRead(x, &xarray));
1077   PetscCall(VecGetArray(y, &yarray));
1078   PetscCall(VecPlaceArray(mpjac->xsub, xarray));
1079   PetscCall(VecPlaceArray(mpjac->ysub, yarray));
1080 
1081   /* apply preconditioner on each matrix block */
1082   PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[0], mpjac->xsub, mpjac->ysub, 0));
1083   PetscCall(KSPSolve(jac->ksp[0], mpjac->xsub, mpjac->ysub));
1084   PetscCall(KSPCheckSolve(jac->ksp[0], pc, mpjac->ysub));
1085   PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[0], mpjac->xsub, mpjac->ysub, 0));
1086   PetscCall(KSPGetConvergedReason(jac->ksp[0], &reason));
1087   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
1088 
1089   PetscCall(VecResetArray(mpjac->xsub));
1090   PetscCall(VecResetArray(mpjac->ysub));
1091   PetscCall(VecRestoreArrayRead(x, &xarray));
1092   PetscCall(VecRestoreArray(y, &yarray));
1093   PetscFunctionReturn(0);
1094 }
1095 
1096 static PetscErrorCode PCMatApply_BJacobi_Multiproc(PC pc, Mat X, Mat Y) {
1097   PC_BJacobi        *jac = (PC_BJacobi *)pc->data;
1098   KSPConvergedReason reason;
1099   Mat                sX, sY;
1100   const PetscScalar *x;
1101   PetscScalar       *y;
1102   PetscInt           m, N, lda, ldb;
1103 
1104   PetscFunctionBegin;
1105   /* apply preconditioner on each matrix block */
1106   PetscCall(MatGetLocalSize(X, &m, NULL));
1107   PetscCall(MatGetSize(X, NULL, &N));
1108   PetscCall(MatDenseGetLDA(X, &lda));
1109   PetscCall(MatDenseGetLDA(Y, &ldb));
1110   PetscCall(MatDenseGetArrayRead(X, &x));
1111   PetscCall(MatDenseGetArrayWrite(Y, &y));
1112   PetscCall(MatCreateDense(PetscObjectComm((PetscObject)jac->ksp[0]), m, PETSC_DECIDE, PETSC_DECIDE, N, (PetscScalar *)x, &sX));
1113   PetscCall(MatCreateDense(PetscObjectComm((PetscObject)jac->ksp[0]), m, PETSC_DECIDE, PETSC_DECIDE, N, y, &sY));
1114   PetscCall(MatDenseSetLDA(sX, lda));
1115   PetscCall(MatDenseSetLDA(sY, ldb));
1116   PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[0], X, Y, 0));
1117   PetscCall(KSPMatSolve(jac->ksp[0], sX, sY));
1118   PetscCall(KSPCheckSolve(jac->ksp[0], pc, NULL));
1119   PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[0], X, Y, 0));
1120   PetscCall(MatDestroy(&sY));
1121   PetscCall(MatDestroy(&sX));
1122   PetscCall(MatDenseRestoreArrayWrite(Y, &y));
1123   PetscCall(MatDenseRestoreArrayRead(X, &x));
1124   PetscCall(KSPGetConvergedReason(jac->ksp[0], &reason));
1125   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
1126   PetscFunctionReturn(0);
1127 }
1128 
1129 static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC pc) {
1130   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
1131   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
1132   PetscInt              m, n;
1133   MPI_Comm              comm, subcomm = 0;
1134   const char           *prefix;
1135   PetscBool             wasSetup = PETSC_TRUE;
1136   VecType               vectype;
1137 
1138   PetscFunctionBegin;
1139   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
1140   PetscCheck(jac->n_local <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Only a single block in a subcommunicator is supported");
1141   jac->n_local = 1; /* currently only a single block is supported for a subcommunicator */
1142   if (!pc->setupcalled) {
1143     wasSetup = PETSC_FALSE;
1144     PetscCall(PetscNewLog(pc, &mpjac));
1145     jac->data = (void *)mpjac;
1146 
1147     /* initialize datastructure mpjac */
1148     if (!jac->psubcomm) {
1149       /* Create default contiguous subcommunicatiors if user does not provide them */
1150       PetscCall(PetscSubcommCreate(comm, &jac->psubcomm));
1151       PetscCall(PetscSubcommSetNumber(jac->psubcomm, jac->n));
1152       PetscCall(PetscSubcommSetType(jac->psubcomm, PETSC_SUBCOMM_CONTIGUOUS));
1153       PetscCall(PetscLogObjectMemory((PetscObject)pc, sizeof(PetscSubcomm)));
1154     }
1155     mpjac->psubcomm = jac->psubcomm;
1156     subcomm         = PetscSubcommChild(mpjac->psubcomm);
1157 
1158     /* Get matrix blocks of pmat */
1159     PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_INITIAL_MATRIX, &mpjac->submats));
1160 
1161     /* create a new PC that processors in each subcomm have copy of */
1162     PetscCall(PetscMalloc1(1, &jac->ksp));
1163     PetscCall(KSPCreate(subcomm, &jac->ksp[0]));
1164     PetscCall(KSPSetErrorIfNotConverged(jac->ksp[0], pc->erroriffailure));
1165     PetscCall(PetscObjectIncrementTabLevel((PetscObject)jac->ksp[0], (PetscObject)pc, 1));
1166     PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)jac->ksp[0]));
1167     PetscCall(KSPSetOperators(jac->ksp[0], mpjac->submats, mpjac->submats));
1168     PetscCall(KSPGetPC(jac->ksp[0], &mpjac->pc));
1169 
1170     PetscCall(PCGetOptionsPrefix(pc, &prefix));
1171     PetscCall(KSPSetOptionsPrefix(jac->ksp[0], prefix));
1172     PetscCall(KSPAppendOptionsPrefix(jac->ksp[0], "sub_"));
1173     PetscCall(KSPGetOptionsPrefix(jac->ksp[0], &prefix));
1174     PetscCall(MatSetOptionsPrefix(mpjac->submats, prefix));
1175 
1176     /* create dummy vectors xsub and ysub */
1177     PetscCall(MatGetLocalSize(mpjac->submats, &m, &n));
1178     PetscCall(VecCreateMPIWithArray(subcomm, 1, n, PETSC_DECIDE, NULL, &mpjac->xsub));
1179     PetscCall(VecCreateMPIWithArray(subcomm, 1, m, PETSC_DECIDE, NULL, &mpjac->ysub));
1180     PetscCall(MatGetVecType(mpjac->submats, &vectype));
1181     PetscCall(VecSetType(mpjac->xsub, vectype));
1182     PetscCall(VecSetType(mpjac->ysub, vectype));
1183     PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)mpjac->xsub));
1184     PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)mpjac->ysub));
1185 
1186     pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiproc;
1187     pc->ops->reset         = PCReset_BJacobi_Multiproc;
1188     pc->ops->destroy       = PCDestroy_BJacobi_Multiproc;
1189     pc->ops->apply         = PCApply_BJacobi_Multiproc;
1190     pc->ops->matapply      = PCMatApply_BJacobi_Multiproc;
1191   } else { /* pc->setupcalled */
1192     subcomm = PetscSubcommChild(mpjac->psubcomm);
1193     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
1194       /* destroy old matrix blocks, then get new matrix blocks */
1195       if (mpjac->submats) PetscCall(MatDestroy(&mpjac->submats));
1196       PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_INITIAL_MATRIX, &mpjac->submats));
1197     } else {
1198       PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_REUSE_MATRIX, &mpjac->submats));
1199     }
1200     PetscCall(KSPSetOperators(jac->ksp[0], mpjac->submats, mpjac->submats));
1201   }
1202 
1203   if (!wasSetup && pc->setfromoptionscalled) PetscCall(KSPSetFromOptions(jac->ksp[0]));
1204   PetscFunctionReturn(0);
1205 }
1206