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