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