xref: /petsc/src/ksp/pc/impls/bjacobi/bjacobi.c (revision 7687e1ecc0da34deafbca176f4acbcceb4e2b94c)
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);CHKERRMPI(ierr);
25   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRMPI(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         = PetscArraycpy(jac->l_lens,jac->g_lens,jac->n_local);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         = PetscArraycpy(jac->l_lens,jac->g_lens+i_start,jac->n_local);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);CHKERRMPI(ierr);
202     if (jac->same_local_solves) {
203       ierr = PetscViewerASCIIPrintf(viewer,"  Local solver is the same for all blocks, as in the following KSP and PC objects on rank 0:\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 = PetscViewerFlush(sviewer);CHKERRQ(ierr);
212         ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
213         ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
214         /*  extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
215         ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
216       } else if (mpjac && jac->ksp && mpjac->psubcomm) {
217         ierr = PetscViewerGetSubViewer(viewer,mpjac->psubcomm->child,&sviewer);CHKERRQ(ierr);
218         if (!mpjac->psubcomm->color) {
219           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
220           ierr = KSPView(*(jac->ksp),sviewer);CHKERRQ(ierr);
221           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
222         }
223         ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr);
224         ierr = PetscViewerRestoreSubViewer(viewer,mpjac->psubcomm->child,&sviewer);CHKERRQ(ierr);
225         ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
226         /*  extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
227         ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
228       } else {
229         ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
230       }
231     } else {
232       PetscInt n_global;
233       ierr = MPIU_Allreduce(&jac->n_local,&n_global,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
234       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
235       ierr = PetscViewerASCIIPrintf(viewer,"  Local solve info for each block is in the following KSP and PC objects:\n");CHKERRQ(ierr);
236       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] number of local blocks = %D, first local block number = %D\n",
237                                                 rank,jac->n_local,jac->first_local);CHKERRQ(ierr);
238       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
239       ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
240       for (i=0; i<jac->n_local; i++) {
241         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] local block number %D\n",rank,i);CHKERRQ(ierr);
242         ierr = KSPView(jac->ksp[i],sviewer);CHKERRQ(ierr);
243         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
244       }
245       ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
246       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
247       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
248       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
249     }
250   } else if (isstring) {
251     ierr = PetscViewerStringSPrintf(viewer," blks=%D",jac->n);CHKERRQ(ierr);
252     ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
253     if (jac->ksp) {ierr = KSPView(jac->ksp[0],sviewer);CHKERRQ(ierr);}
254     ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
255   } else if (isdraw) {
256     PetscDraw draw;
257     char      str[25];
258     PetscReal x,y,bottom,h;
259 
260     ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
261     ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
262     ierr   = PetscSNPrintf(str,25,"Number blocks %D",jac->n);CHKERRQ(ierr);
263     ierr   = PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr);
264     bottom = y - h;
265     ierr   = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
266     /* warning the communicator on viewer is different then on ksp in parallel */
267     if (jac->ksp) {ierr = KSPView(jac->ksp[0],viewer);CHKERRQ(ierr);}
268     ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
269   }
270   PetscFunctionReturn(0);
271 }
272 
273 /* -------------------------------------------------------------------------------------*/
274 
275 static PetscErrorCode  PCBJacobiGetSubKSP_BJacobi(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp)
276 {
277   PC_BJacobi *jac = (PC_BJacobi*)pc->data;
278 
279   PetscFunctionBegin;
280   if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call KSPSetUp() or PCSetUp() first");
281 
282   if (n_local) *n_local = jac->n_local;
283   if (first_local) *first_local = jac->first_local;
284   *ksp                   = jac->ksp;
285   jac->same_local_solves = PETSC_FALSE;        /* Assume that local solves are now different;
286                                                   not necessarily true though!  This flag is
287                                                   used only for PCView_BJacobi() */
288   PetscFunctionReturn(0);
289 }
290 
291 static PetscErrorCode  PCBJacobiSetTotalBlocks_BJacobi(PC pc,PetscInt blocks,PetscInt *lens)
292 {
293   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;
294   PetscErrorCode ierr;
295 
296   PetscFunctionBegin;
297   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");
298   jac->n = blocks;
299   if (!lens) jac->g_lens = NULL;
300   else {
301     ierr = PetscMalloc1(blocks,&jac->g_lens);CHKERRQ(ierr);
302     ierr = PetscLogObjectMemory((PetscObject)pc,blocks*sizeof(PetscInt));CHKERRQ(ierr);
303     ierr = PetscArraycpy(jac->g_lens,lens,blocks);CHKERRQ(ierr);
304   }
305   PetscFunctionReturn(0);
306 }
307 
308 static PetscErrorCode  PCBJacobiGetTotalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[])
309 {
310   PC_BJacobi *jac = (PC_BJacobi*) pc->data;
311 
312   PetscFunctionBegin;
313   *blocks = jac->n;
314   if (lens) *lens = jac->g_lens;
315   PetscFunctionReturn(0);
316 }
317 
318 static PetscErrorCode  PCBJacobiSetLocalBlocks_BJacobi(PC pc,PetscInt blocks,const PetscInt lens[])
319 {
320   PC_BJacobi     *jac;
321   PetscErrorCode ierr;
322 
323   PetscFunctionBegin;
324   jac = (PC_BJacobi*)pc->data;
325 
326   jac->n_local = blocks;
327   if (!lens) jac->l_lens = NULL;
328   else {
329     ierr = PetscMalloc1(blocks,&jac->l_lens);CHKERRQ(ierr);
330     ierr = PetscLogObjectMemory((PetscObject)pc,blocks*sizeof(PetscInt));CHKERRQ(ierr);
331     ierr = PetscArraycpy(jac->l_lens,lens,blocks);CHKERRQ(ierr);
332   }
333   PetscFunctionReturn(0);
334 }
335 
336 static PetscErrorCode  PCBJacobiGetLocalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[])
337 {
338   PC_BJacobi *jac = (PC_BJacobi*) pc->data;
339 
340   PetscFunctionBegin;
341   *blocks = jac->n_local;
342   if (lens) *lens = jac->l_lens;
343   PetscFunctionReturn(0);
344 }
345 
346 /* -------------------------------------------------------------------------------------*/
347 
348 /*@C
349    PCBJacobiGetSubKSP - Gets the local KSP contexts for all blocks on
350    this processor.
351 
352    Not Collective
353 
354    Input Parameter:
355 .  pc - the preconditioner context
356 
357    Output Parameters:
358 +  n_local - the number of blocks on this processor, or NULL
359 .  first_local - the global number of the first block on this processor, or NULL
360 -  ksp - the array of KSP contexts
361 
362    Notes:
363    After PCBJacobiGetSubKSP() the array of KSP contexts is not to be freed.
364 
365    Currently for some matrix implementations only 1 block per processor
366    is supported.
367 
368    You must call KSPSetUp() or PCSetUp() before calling PCBJacobiGetSubKSP().
369 
370    Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs.
371       You can call PCBJacobiGetSubKSP(pc,nlocal,firstlocal,PETSC_NULL_KSP,ierr) to determine how large the
372       KSP array must be.
373 
374    Level: advanced
375 
376 .seealso: PCASMGetSubKSP()
377 @*/
378 PetscErrorCode  PCBJacobiGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
379 {
380   PetscErrorCode ierr;
381 
382   PetscFunctionBegin;
383   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
384   ierr = PetscUseMethod(pc,"PCBJacobiGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr);
385   PetscFunctionReturn(0);
386 }
387 
388 /*@
389    PCBJacobiSetTotalBlocks - Sets the global number of blocks for the block
390    Jacobi preconditioner.
391 
392    Collective on PC
393 
394    Input Parameters:
395 +  pc - the preconditioner context
396 .  blocks - the number of blocks
397 -  lens - [optional] integer array containing the size of each block
398 
399    Options Database Key:
400 .  -pc_bjacobi_blocks <blocks> - Sets the number of global blocks
401 
402    Notes:
403    Currently only a limited number of blocking configurations are supported.
404    All processors sharing the PC must call this routine with the same data.
405 
406    Level: intermediate
407 
408 .seealso: PCSetUseAmat(), PCBJacobiSetLocalBlocks()
409 @*/
410 PetscErrorCode  PCBJacobiSetTotalBlocks(PC pc,PetscInt blocks,const PetscInt lens[])
411 {
412   PetscErrorCode ierr;
413 
414   PetscFunctionBegin;
415   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
416   if (blocks <= 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Must have positive blocks");
417   ierr = PetscTryMethod(pc,"PCBJacobiSetTotalBlocks_C",(PC,PetscInt,const PetscInt[]),(pc,blocks,lens));CHKERRQ(ierr);
418   PetscFunctionReturn(0);
419 }
420 
421 /*@C
422    PCBJacobiGetTotalBlocks - Gets the global number of blocks for the block
423    Jacobi preconditioner.
424 
425    Not Collective
426 
427    Input Parameter:
428 .  pc - the preconditioner context
429 
430    Output parameters:
431 +  blocks - the number of blocks
432 -  lens - integer array containing the size of each block
433 
434    Level: intermediate
435 
436 .seealso: PCSetUseAmat(), PCBJacobiGetLocalBlocks()
437 @*/
438 PetscErrorCode  PCBJacobiGetTotalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
439 {
440   PetscErrorCode ierr;
441 
442   PetscFunctionBegin;
443   PetscValidHeaderSpecific(pc, PC_CLASSID,1);
444   PetscValidIntPointer(blocks,2);
445   ierr = PetscUseMethod(pc,"PCBJacobiGetTotalBlocks_C",(PC,PetscInt*, const PetscInt *[]),(pc,blocks,lens));CHKERRQ(ierr);
446   PetscFunctionReturn(0);
447 }
448 
449 /*@
450    PCBJacobiSetLocalBlocks - Sets the local number of blocks for the block
451    Jacobi preconditioner.
452 
453    Not Collective
454 
455    Input Parameters:
456 +  pc - the preconditioner context
457 .  blocks - the number of blocks
458 -  lens - [optional] integer array containing size of each block
459 
460    Options Database Key:
461 .  -pc_bjacobi_local_blocks <blocks> - Sets the number of local blocks
462 
463    Note:
464    Currently only a limited number of blocking configurations are supported.
465 
466    Level: intermediate
467 
468 .seealso: PCSetUseAmat(), PCBJacobiSetTotalBlocks()
469 @*/
470 PetscErrorCode  PCBJacobiSetLocalBlocks(PC pc,PetscInt blocks,const PetscInt lens[])
471 {
472   PetscErrorCode ierr;
473 
474   PetscFunctionBegin;
475   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
476   if (blocks < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have nonegative blocks");
477   ierr = PetscTryMethod(pc,"PCBJacobiSetLocalBlocks_C",(PC,PetscInt,const PetscInt []),(pc,blocks,lens));CHKERRQ(ierr);
478   PetscFunctionReturn(0);
479 }
480 
481 /*@C
482    PCBJacobiGetLocalBlocks - Gets the local number of blocks for the block
483    Jacobi preconditioner.
484 
485    Not Collective
486 
487    Input Parameters:
488 +  pc - the preconditioner context
489 .  blocks - the number of blocks
490 -  lens - [optional] integer array containing size of each block
491 
492    Note:
493    Currently only a limited number of blocking configurations are supported.
494 
495    Level: intermediate
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 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
541            PCASM, PCSetUseAmat(), PCGetUseAmat(), PCBJacobiGetSubKSP(), PCBJacobiSetTotalBlocks(),
542            PCBJacobiSetLocalBlocks(), PCSetModifySubMatrices()
543 M*/
544 
545 PETSC_EXTERN PetscErrorCode PCCreate_BJacobi(PC pc)
546 {
547   PetscErrorCode ierr;
548   PetscMPIInt    rank;
549   PC_BJacobi     *jac;
550 
551   PetscFunctionBegin;
552   ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr);
553   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRMPI(ierr);
554 
555   pc->ops->apply           = NULL;
556   pc->ops->matapply        = NULL;
557   pc->ops->applytranspose  = NULL;
558   pc->ops->setup           = PCSetUp_BJacobi;
559   pc->ops->destroy         = PCDestroy_BJacobi;
560   pc->ops->setfromoptions  = PCSetFromOptions_BJacobi;
561   pc->ops->view            = PCView_BJacobi;
562   pc->ops->applyrichardson = NULL;
563 
564   pc->data               = (void*)jac;
565   jac->n                 = -1;
566   jac->n_local           = -1;
567   jac->first_local       = rank;
568   jac->ksp               = NULL;
569   jac->same_local_solves = PETSC_TRUE;
570   jac->g_lens            = NULL;
571   jac->l_lens            = NULL;
572   jac->psubcomm          = NULL;
573 
574   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiGetSubKSP_C",PCBJacobiGetSubKSP_BJacobi);CHKERRQ(ierr);
575   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiSetTotalBlocks_C",PCBJacobiSetTotalBlocks_BJacobi);CHKERRQ(ierr);
576   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiGetTotalBlocks_C",PCBJacobiGetTotalBlocks_BJacobi);CHKERRQ(ierr);
577   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiSetLocalBlocks_C",PCBJacobiSetLocalBlocks_BJacobi);CHKERRQ(ierr);
578   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiGetLocalBlocks_C",PCBJacobiGetLocalBlocks_BJacobi);CHKERRQ(ierr);
579   PetscFunctionReturn(0);
580 }
581 
582 /* --------------------------------------------------------------------------------------------*/
583 /*
584         These are for a single block per processor; works for AIJ, BAIJ; Seq and MPI
585 */
586 static PetscErrorCode PCReset_BJacobi_Singleblock(PC pc)
587 {
588   PC_BJacobi             *jac  = (PC_BJacobi*)pc->data;
589   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
590   PetscErrorCode         ierr;
591 
592   PetscFunctionBegin;
593   ierr = KSPReset(jac->ksp[0]);CHKERRQ(ierr);
594   ierr = VecDestroy(&bjac->x);CHKERRQ(ierr);
595   ierr = VecDestroy(&bjac->y);CHKERRQ(ierr);
596   PetscFunctionReturn(0);
597 }
598 
599 static PetscErrorCode PCDestroy_BJacobi_Singleblock(PC pc)
600 {
601   PC_BJacobi             *jac  = (PC_BJacobi*)pc->data;
602   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
603   PetscErrorCode         ierr;
604 
605   PetscFunctionBegin;
606   ierr = PCReset_BJacobi_Singleblock(pc);CHKERRQ(ierr);
607   ierr = KSPDestroy(&jac->ksp[0]);CHKERRQ(ierr);
608   ierr = PetscFree(jac->ksp);CHKERRQ(ierr);
609   ierr = PetscFree(jac->l_lens);CHKERRQ(ierr);
610   ierr = PetscFree(jac->g_lens);CHKERRQ(ierr);
611   ierr = PetscFree(bjac);CHKERRQ(ierr);
612   ierr = PetscFree(pc->data);CHKERRQ(ierr);
613   PetscFunctionReturn(0);
614 }
615 
616 static PetscErrorCode PCSetUpOnBlocks_BJacobi_Singleblock(PC pc)
617 {
618   PetscErrorCode     ierr;
619   PC_BJacobi         *jac = (PC_BJacobi*)pc->data;
620   KSP                subksp = jac->ksp[0];
621   KSPConvergedReason reason;
622 
623   PetscFunctionBegin;
624   ierr = KSPSetUp(subksp);CHKERRQ(ierr);
625   ierr = KSPGetConvergedReason(subksp,&reason);CHKERRQ(ierr);
626   if (reason == KSP_DIVERGED_PC_FAILED) {
627     pc->failedreason = PC_SUBPC_ERROR;
628   }
629   PetscFunctionReturn(0);
630 }
631 
632 static PetscErrorCode PCApply_BJacobi_Singleblock(PC pc,Vec x,Vec y)
633 {
634   PetscErrorCode         ierr;
635   PC_BJacobi             *jac  = (PC_BJacobi*)pc->data;
636   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
637 
638   PetscFunctionBegin;
639   ierr = VecGetLocalVectorRead(x, bjac->x);CHKERRQ(ierr);
640   ierr = VecGetLocalVector(y, bjac->y);CHKERRQ(ierr);
641  /* Since the inner KSP matrix may point directly to the diagonal block of an MPI matrix the inner
642      matrix may change even if the outer KSP/PC has not updated the preconditioner, this will trigger a rebuild
643      of the inner preconditioner automatically unless we pass down the outer preconditioners reuse flag.*/
644   ierr = KSPSetReusePreconditioner(jac->ksp[0],pc->reusepreconditioner);CHKERRQ(ierr);
645   ierr = KSPSolve(jac->ksp[0],bjac->x,bjac->y);CHKERRQ(ierr);
646   ierr = KSPCheckSolve(jac->ksp[0],pc,bjac->y);CHKERRQ(ierr);
647   ierr = VecRestoreLocalVectorRead(x, bjac->x);CHKERRQ(ierr);
648   ierr = VecRestoreLocalVector(y, bjac->y);CHKERRQ(ierr);
649   PetscFunctionReturn(0);
650 }
651 
652 static PetscErrorCode PCMatApply_BJacobi_Singleblock(PC pc,Mat X,Mat Y)
653 {
654   PC_BJacobi     *jac  = (PC_BJacobi*)pc->data;
655   Mat            sX,sY;
656   PetscErrorCode ierr;
657 
658   PetscFunctionBegin;
659  /* Since the inner KSP matrix may point directly to the diagonal block of an MPI matrix the inner
660      matrix may change even if the outer KSP/PC has not updated the preconditioner, this will trigger a rebuild
661      of the inner preconditioner automatically unless we pass down the outer preconditioners reuse flag.*/
662   ierr = KSPSetReusePreconditioner(jac->ksp[0],pc->reusepreconditioner);CHKERRQ(ierr);
663   ierr = MatDenseGetLocalMatrix(X,&sX);CHKERRQ(ierr);
664   ierr = MatDenseGetLocalMatrix(Y,&sY);CHKERRQ(ierr);
665   ierr = KSPMatSolve(jac->ksp[0],sX,sY);CHKERRQ(ierr);
666   PetscFunctionReturn(0);
667 }
668 
669 static PetscErrorCode PCApplySymmetricLeft_BJacobi_Singleblock(PC pc,Vec x,Vec y)
670 {
671   PetscErrorCode         ierr;
672   PC_BJacobi             *jac  = (PC_BJacobi*)pc->data;
673   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
674   PetscScalar            *y_array;
675   const PetscScalar      *x_array;
676   PC                     subpc;
677 
678   PetscFunctionBegin;
679   /*
680       The VecPlaceArray() is to avoid having to copy the
681     y vector into the bjac->x vector. The reason for
682     the bjac->x vector is that we need a sequential vector
683     for the sequential solve.
684   */
685   ierr = VecGetArrayRead(x,&x_array);CHKERRQ(ierr);
686   ierr = VecGetArray(y,&y_array);CHKERRQ(ierr);
687   ierr = VecPlaceArray(bjac->x,x_array);CHKERRQ(ierr);
688   ierr = VecPlaceArray(bjac->y,y_array);CHKERRQ(ierr);
689   /* apply the symmetric left portion of the inner PC operator */
690   /* note this by-passes the inner KSP and its options completely */
691   ierr = KSPGetPC(jac->ksp[0],&subpc);CHKERRQ(ierr);
692   ierr = PCApplySymmetricLeft(subpc,bjac->x,bjac->y);CHKERRQ(ierr);
693   ierr = VecResetArray(bjac->x);CHKERRQ(ierr);
694   ierr = VecResetArray(bjac->y);CHKERRQ(ierr);
695   ierr = VecRestoreArrayRead(x,&x_array);CHKERRQ(ierr);
696   ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr);
697   PetscFunctionReturn(0);
698 }
699 
700 static PetscErrorCode PCApplySymmetricRight_BJacobi_Singleblock(PC pc,Vec x,Vec y)
701 {
702   PetscErrorCode         ierr;
703   PC_BJacobi             *jac  = (PC_BJacobi*)pc->data;
704   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
705   PetscScalar            *y_array;
706   const PetscScalar      *x_array;
707   PC                     subpc;
708 
709   PetscFunctionBegin;
710   /*
711       The VecPlaceArray() is to avoid having to copy the
712     y vector into the bjac->x vector. The reason for
713     the bjac->x vector is that we need a sequential vector
714     for the sequential solve.
715   */
716   ierr = VecGetArrayRead(x,&x_array);CHKERRQ(ierr);
717   ierr = VecGetArray(y,&y_array);CHKERRQ(ierr);
718   ierr = VecPlaceArray(bjac->x,x_array);CHKERRQ(ierr);
719   ierr = VecPlaceArray(bjac->y,y_array);CHKERRQ(ierr);
720 
721   /* apply the symmetric right portion of the inner PC operator */
722   /* note this by-passes the inner KSP and its options completely */
723 
724   ierr = KSPGetPC(jac->ksp[0],&subpc);CHKERRQ(ierr);
725   ierr = PCApplySymmetricRight(subpc,bjac->x,bjac->y);CHKERRQ(ierr);
726 
727   ierr = VecRestoreArrayRead(x,&x_array);CHKERRQ(ierr);
728   ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr);
729   PetscFunctionReturn(0);
730 }
731 
732 static PetscErrorCode PCApplyTranspose_BJacobi_Singleblock(PC pc,Vec x,Vec y)
733 {
734   PetscErrorCode         ierr;
735   PC_BJacobi             *jac  = (PC_BJacobi*)pc->data;
736   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
737   PetscScalar            *y_array;
738   const PetscScalar      *x_array;
739 
740   PetscFunctionBegin;
741   /*
742       The VecPlaceArray() is to avoid having to copy the
743     y vector into the bjac->x vector. The reason for
744     the bjac->x vector is that we need a sequential vector
745     for the sequential solve.
746   */
747   ierr = VecGetArrayRead(x,&x_array);CHKERRQ(ierr);
748   ierr = VecGetArray(y,&y_array);CHKERRQ(ierr);
749   ierr = VecPlaceArray(bjac->x,x_array);CHKERRQ(ierr);
750   ierr = VecPlaceArray(bjac->y,y_array);CHKERRQ(ierr);
751   ierr = KSPSolveTranspose(jac->ksp[0],bjac->x,bjac->y);CHKERRQ(ierr);
752   ierr = KSPCheckSolve(jac->ksp[0],pc,bjac->y);CHKERRQ(ierr);
753   ierr = VecResetArray(bjac->x);CHKERRQ(ierr);
754   ierr = VecResetArray(bjac->y);CHKERRQ(ierr);
755   ierr = VecRestoreArrayRead(x,&x_array);CHKERRQ(ierr);
756   ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr);
757   PetscFunctionReturn(0);
758 }
759 
760 static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC pc,Mat mat,Mat pmat)
761 {
762   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
763   PetscErrorCode         ierr;
764   PetscInt               m;
765   KSP                    ksp;
766   PC_BJacobi_Singleblock *bjac;
767   PetscBool              wasSetup = PETSC_TRUE;
768   VecType                vectype;
769 
770   PetscFunctionBegin;
771   if (!pc->setupcalled) {
772     const char *prefix;
773 
774     if (!jac->ksp) {
775       wasSetup = PETSC_FALSE;
776 
777       ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr);
778       ierr = KSPSetErrorIfNotConverged(ksp,pc->erroriffailure);CHKERRQ(ierr);
779       ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
780       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);CHKERRQ(ierr);
781       ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
782       ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
783       ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
784       ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
785 
786       pc->ops->reset               = PCReset_BJacobi_Singleblock;
787       pc->ops->destroy             = PCDestroy_BJacobi_Singleblock;
788       pc->ops->apply               = PCApply_BJacobi_Singleblock;
789       pc->ops->matapply            = PCMatApply_BJacobi_Singleblock;
790       pc->ops->applysymmetricleft  = PCApplySymmetricLeft_BJacobi_Singleblock;
791       pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Singleblock;
792       pc->ops->applytranspose      = PCApplyTranspose_BJacobi_Singleblock;
793       pc->ops->setuponblocks       = PCSetUpOnBlocks_BJacobi_Singleblock;
794 
795       ierr        = PetscMalloc1(1,&jac->ksp);CHKERRQ(ierr);
796       jac->ksp[0] = ksp;
797 
798       ierr      = PetscNewLog(pc,&bjac);CHKERRQ(ierr);
799       jac->data = (void*)bjac;
800     } else {
801       ksp  = jac->ksp[0];
802       bjac = (PC_BJacobi_Singleblock*)jac->data;
803     }
804 
805     /*
806       The reason we need to generate these vectors is to serve
807       as the right-hand side and solution vector for the solve on the
808       block. We do not need to allocate space for the vectors since
809       that is provided via VecPlaceArray() just before the call to
810       KSPSolve() on the block.
811     */
812     ierr = MatGetSize(pmat,&m,&m);CHKERRQ(ierr);
813     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,m,NULL,&bjac->x);CHKERRQ(ierr);
814     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,m,NULL,&bjac->y);CHKERRQ(ierr);
815     ierr = MatGetVecType(pmat,&vectype);CHKERRQ(ierr);
816     ierr = VecSetType(bjac->x,vectype);CHKERRQ(ierr);
817     ierr = VecSetType(bjac->y,vectype);CHKERRQ(ierr);
818     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)bjac->x);CHKERRQ(ierr);
819     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)bjac->y);CHKERRQ(ierr);
820   } else {
821     ksp  = jac->ksp[0];
822     bjac = (PC_BJacobi_Singleblock*)jac->data;
823   }
824   if (pc->useAmat) {
825     ierr = KSPSetOperators(ksp,mat,pmat);CHKERRQ(ierr);
826   } else {
827     ierr = KSPSetOperators(ksp,pmat,pmat);CHKERRQ(ierr);
828   }
829   if (!wasSetup && pc->setfromoptionscalled) {
830     /* If PCSetFromOptions_BJacobi is called later, KSPSetFromOptions will be called at that time. */
831     ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
832   }
833   PetscFunctionReturn(0);
834 }
835 
836 /* ---------------------------------------------------------------------------------------------*/
837 static PetscErrorCode PCReset_BJacobi_Multiblock(PC pc)
838 {
839   PC_BJacobi            *jac  = (PC_BJacobi*)pc->data;
840   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
841   PetscErrorCode        ierr;
842   PetscInt              i;
843 
844   PetscFunctionBegin;
845   if (bjac && bjac->pmat) {
846     ierr = MatDestroyMatrices(jac->n_local,&bjac->pmat);CHKERRQ(ierr);
847     if (pc->useAmat) {
848       ierr = MatDestroyMatrices(jac->n_local,&bjac->mat);CHKERRQ(ierr);
849     }
850   }
851 
852   for (i=0; i<jac->n_local; i++) {
853     ierr = KSPReset(jac->ksp[i]);CHKERRQ(ierr);
854     if (bjac && bjac->x) {
855       ierr = VecDestroy(&bjac->x[i]);CHKERRQ(ierr);
856       ierr = VecDestroy(&bjac->y[i]);CHKERRQ(ierr);
857       ierr = ISDestroy(&bjac->is[i]);CHKERRQ(ierr);
858     }
859   }
860   ierr = PetscFree(jac->l_lens);CHKERRQ(ierr);
861   ierr = PetscFree(jac->g_lens);CHKERRQ(ierr);
862   PetscFunctionReturn(0);
863 }
864 
865 static PetscErrorCode PCDestroy_BJacobi_Multiblock(PC pc)
866 {
867   PC_BJacobi            *jac  = (PC_BJacobi*)pc->data;
868   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
869   PetscErrorCode        ierr;
870   PetscInt              i;
871 
872   PetscFunctionBegin;
873   ierr = PCReset_BJacobi_Multiblock(pc);CHKERRQ(ierr);
874   if (bjac) {
875     ierr = PetscFree2(bjac->x,bjac->y);CHKERRQ(ierr);
876     ierr = PetscFree(bjac->starts);CHKERRQ(ierr);
877     ierr = PetscFree(bjac->is);CHKERRQ(ierr);
878   }
879   ierr = PetscFree(jac->data);CHKERRQ(ierr);
880   for (i=0; i<jac->n_local; i++) {
881     ierr = KSPDestroy(&jac->ksp[i]);CHKERRQ(ierr);
882   }
883   ierr = PetscFree(jac->ksp);CHKERRQ(ierr);
884   ierr = PetscFree(pc->data);CHKERRQ(ierr);
885   PetscFunctionReturn(0);
886 }
887 
888 static PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiblock(PC pc)
889 {
890   PC_BJacobi         *jac = (PC_BJacobi*)pc->data;
891   PetscErrorCode     ierr;
892   PetscInt           i,n_local = jac->n_local;
893   KSPConvergedReason reason;
894 
895   PetscFunctionBegin;
896   for (i=0; i<n_local; i++) {
897     ierr = KSPSetUp(jac->ksp[i]);CHKERRQ(ierr);
898     ierr = KSPGetConvergedReason(jac->ksp[i],&reason);CHKERRQ(ierr);
899     if (reason == KSP_DIVERGED_PC_FAILED) {
900       pc->failedreason = PC_SUBPC_ERROR;
901     }
902   }
903   PetscFunctionReturn(0);
904 }
905 
906 /*
907       Preconditioner for block Jacobi
908 */
909 static PetscErrorCode PCApply_BJacobi_Multiblock(PC pc,Vec x,Vec y)
910 {
911   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
912   PetscErrorCode        ierr;
913   PetscInt              i,n_local = jac->n_local;
914   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
915   PetscScalar           *yin;
916   const PetscScalar     *xin;
917 
918   PetscFunctionBegin;
919   ierr = VecGetArrayRead(x,&xin);CHKERRQ(ierr);
920   ierr = VecGetArray(y,&yin);CHKERRQ(ierr);
921   for (i=0; i<n_local; i++) {
922     /*
923        To avoid copying the subvector from x into a workspace we instead
924        make the workspace vector array point to the subpart of the array of
925        the global vector.
926     */
927     ierr = VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);CHKERRQ(ierr);
928     ierr = VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);CHKERRQ(ierr);
929 
930     ierr = PetscLogEventBegin(PC_ApplyOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);CHKERRQ(ierr);
931     ierr = KSPSolve(jac->ksp[i],bjac->x[i],bjac->y[i]);CHKERRQ(ierr);
932     ierr = KSPCheckSolve(jac->ksp[i],pc,bjac->y[i]);CHKERRQ(ierr);
933     ierr = PetscLogEventEnd(PC_ApplyOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);CHKERRQ(ierr);
934 
935     ierr = VecResetArray(bjac->x[i]);CHKERRQ(ierr);
936     ierr = VecResetArray(bjac->y[i]);CHKERRQ(ierr);
937   }
938   ierr = VecRestoreArrayRead(x,&xin);CHKERRQ(ierr);
939   ierr = VecRestoreArray(y,&yin);CHKERRQ(ierr);
940   PetscFunctionReturn(0);
941 }
942 
943 /*
944       Preconditioner for block Jacobi
945 */
946 static PetscErrorCode PCApplyTranspose_BJacobi_Multiblock(PC pc,Vec x,Vec y)
947 {
948   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
949   PetscErrorCode        ierr;
950   PetscInt              i,n_local = jac->n_local;
951   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
952   PetscScalar           *yin;
953   const PetscScalar     *xin;
954 
955   PetscFunctionBegin;
956   ierr = VecGetArrayRead(x,&xin);CHKERRQ(ierr);
957   ierr = VecGetArray(y,&yin);CHKERRQ(ierr);
958   for (i=0; i<n_local; i++) {
959     /*
960        To avoid copying the subvector from x into a workspace we instead
961        make the workspace vector array point to the subpart of the array of
962        the global vector.
963     */
964     ierr = VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);CHKERRQ(ierr);
965     ierr = VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);CHKERRQ(ierr);
966 
967     ierr = PetscLogEventBegin(PC_ApplyTransposeOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);CHKERRQ(ierr);
968     ierr = KSPSolveTranspose(jac->ksp[i],bjac->x[i],bjac->y[i]);CHKERRQ(ierr);
969     ierr = KSPCheckSolve(jac->ksp[i],pc,bjac->y[i]);CHKERRQ(ierr);
970     ierr = PetscLogEventEnd(PC_ApplyTransposeOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);CHKERRQ(ierr);
971 
972     ierr = VecResetArray(bjac->x[i]);CHKERRQ(ierr);
973     ierr = VecResetArray(bjac->y[i]);CHKERRQ(ierr);
974   }
975   ierr = VecRestoreArrayRead(x,&xin);CHKERRQ(ierr);
976   ierr = VecRestoreArray(y,&yin);CHKERRQ(ierr);
977   PetscFunctionReturn(0);
978 }
979 
980 static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC pc,Mat mat,Mat pmat)
981 {
982   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
983   PetscErrorCode        ierr;
984   PetscInt              m,n_local,N,M,start,i;
985   const char            *prefix,*pprefix,*mprefix;
986   KSP                   ksp;
987   Vec                   x,y;
988   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
989   PC                    subpc;
990   IS                    is;
991   MatReuse              scall;
992   VecType               vectype;
993 
994   PetscFunctionBegin;
995   ierr = MatGetLocalSize(pc->pmat,&M,&N);CHKERRQ(ierr);
996 
997   n_local = jac->n_local;
998 
999   if (pc->useAmat) {
1000     PetscBool same;
1001     ierr = PetscObjectTypeCompare((PetscObject)mat,((PetscObject)pmat)->type_name,&same);CHKERRQ(ierr);
1002     if (!same) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Matrices not of same type");
1003   }
1004 
1005   if (!pc->setupcalled) {
1006     scall = MAT_INITIAL_MATRIX;
1007 
1008     if (!jac->ksp) {
1009       pc->ops->reset         = PCReset_BJacobi_Multiblock;
1010       pc->ops->destroy       = PCDestroy_BJacobi_Multiblock;
1011       pc->ops->apply         = PCApply_BJacobi_Multiblock;
1012       pc->ops->matapply      = NULL;
1013       pc->ops->applytranspose= PCApplyTranspose_BJacobi_Multiblock;
1014       pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiblock;
1015 
1016       ierr = PetscNewLog(pc,&bjac);CHKERRQ(ierr);
1017       ierr = PetscMalloc1(n_local,&jac->ksp);CHKERRQ(ierr);
1018       ierr = PetscLogObjectMemory((PetscObject)pc,sizeof(n_local*sizeof(KSP)));CHKERRQ(ierr);
1019       ierr = PetscMalloc2(n_local,&bjac->x,n_local,&bjac->y);CHKERRQ(ierr);
1020       ierr = PetscMalloc1(n_local,&bjac->starts);CHKERRQ(ierr);
1021       ierr = PetscLogObjectMemory((PetscObject)pc,sizeof(n_local*sizeof(PetscScalar)));CHKERRQ(ierr);
1022 
1023       jac->data = (void*)bjac;
1024       ierr      = PetscMalloc1(n_local,&bjac->is);CHKERRQ(ierr);
1025       ierr      = PetscLogObjectMemory((PetscObject)pc,sizeof(n_local*sizeof(IS)));CHKERRQ(ierr);
1026 
1027       for (i=0; i<n_local; i++) {
1028         ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr);
1029         ierr = KSPSetErrorIfNotConverged(ksp,pc->erroriffailure);CHKERRQ(ierr);
1030         ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1031         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);CHKERRQ(ierr);
1032         ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
1033         ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr);
1034         ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
1035         ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
1036         ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
1037 
1038         jac->ksp[i] = ksp;
1039       }
1040     } else {
1041       bjac = (PC_BJacobi_Multiblock*)jac->data;
1042     }
1043 
1044     start = 0;
1045     ierr = MatGetVecType(pmat,&vectype);CHKERRQ(ierr);
1046     for (i=0; i<n_local; i++) {
1047       m = jac->l_lens[i];
1048       /*
1049       The reason we need to generate these vectors is to serve
1050       as the right-hand side and solution vector for the solve on the
1051       block. We do not need to allocate space for the vectors since
1052       that is provided via VecPlaceArray() just before the call to
1053       KSPSolve() on the block.
1054 
1055       */
1056       ierr = VecCreateSeq(PETSC_COMM_SELF,m,&x);CHKERRQ(ierr);
1057       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,m,NULL,&y);CHKERRQ(ierr);
1058       ierr = VecSetType(x,vectype);CHKERRQ(ierr);
1059       ierr = VecSetType(y,vectype);CHKERRQ(ierr);
1060       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)x);CHKERRQ(ierr);
1061       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)y);CHKERRQ(ierr);
1062 
1063       bjac->x[i]      = x;
1064       bjac->y[i]      = y;
1065       bjac->starts[i] = start;
1066 
1067       ierr        = ISCreateStride(PETSC_COMM_SELF,m,start,1,&is);CHKERRQ(ierr);
1068       bjac->is[i] = is;
1069       ierr        = PetscLogObjectParent((PetscObject)pc,(PetscObject)is);CHKERRQ(ierr);
1070 
1071       start += m;
1072     }
1073   } else {
1074     bjac = (PC_BJacobi_Multiblock*)jac->data;
1075     /*
1076        Destroy the blocks from the previous iteration
1077     */
1078     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
1079       ierr = MatDestroyMatrices(n_local,&bjac->pmat);CHKERRQ(ierr);
1080       if (pc->useAmat) {
1081         ierr = MatDestroyMatrices(n_local,&bjac->mat);CHKERRQ(ierr);
1082       }
1083       scall = MAT_INITIAL_MATRIX;
1084     } else scall = MAT_REUSE_MATRIX;
1085   }
1086 
1087   ierr = MatCreateSubMatrices(pmat,n_local,bjac->is,bjac->is,scall,&bjac->pmat);CHKERRQ(ierr);
1088   if (pc->useAmat) {
1089     ierr = MatCreateSubMatrices(mat,n_local,bjac->is,bjac->is,scall,&bjac->mat);CHKERRQ(ierr);
1090   }
1091   /* Return control to the user so that the submatrices can be modified (e.g., to apply
1092      different boundary conditions for the submatrices than for the global problem) */
1093   ierr = PCModifySubMatrices(pc,n_local,bjac->is,bjac->is,bjac->pmat,pc->modifysubmatricesP);CHKERRQ(ierr);
1094 
1095   ierr = PetscObjectGetOptionsPrefix((PetscObject)pmat,&pprefix);CHKERRQ(ierr);
1096   for (i=0; i<n_local; i++) {
1097     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)bjac->pmat[i]);CHKERRQ(ierr);
1098     ierr = PetscObjectSetOptionsPrefix((PetscObject)bjac->pmat[i],pprefix);CHKERRQ(ierr);
1099     ierr = KSPGetOptionsPrefix(jac->ksp[i],&prefix);CHKERRQ(ierr);
1100     if (pc->useAmat) {
1101       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)bjac->mat[i]);CHKERRQ(ierr);
1102       ierr = PetscObjectGetOptionsPrefix((PetscObject)mat,&mprefix);CHKERRQ(ierr);
1103       ierr = PetscObjectSetOptionsPrefix((PetscObject)bjac->mat[i],mprefix);CHKERRQ(ierr);
1104       ierr = KSPSetOperators(jac->ksp[i],bjac->mat[i],bjac->pmat[i]);CHKERRQ(ierr);
1105       ierr = MatSetOptionsPrefix(bjac->mat[i],prefix);CHKERRQ(ierr);
1106     } else {
1107       ierr = KSPSetOperators(jac->ksp[i],bjac->pmat[i],bjac->pmat[i]);CHKERRQ(ierr);
1108     }
1109     ierr = MatSetOptionsPrefix(bjac->pmat[i],prefix);CHKERRQ(ierr);
1110     if (pc->setfromoptionscalled) {
1111       ierr = KSPSetFromOptions(jac->ksp[i]);CHKERRQ(ierr);
1112     }
1113   }
1114   PetscFunctionReturn(0);
1115 }
1116 
1117 /* ---------------------------------------------------------------------------------------------*/
1118 /*
1119       These are for a single block with multiple processes;
1120 */
1121 static PetscErrorCode PCReset_BJacobi_Multiproc(PC pc)
1122 {
1123   PC_BJacobi           *jac   = (PC_BJacobi*)pc->data;
1124   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1125   PetscErrorCode       ierr;
1126 
1127   PetscFunctionBegin;
1128   ierr = VecDestroy(&mpjac->ysub);CHKERRQ(ierr);
1129   ierr = VecDestroy(&mpjac->xsub);CHKERRQ(ierr);
1130   ierr = MatDestroy(&mpjac->submats);CHKERRQ(ierr);
1131   if (jac->ksp) {ierr = KSPReset(jac->ksp[0]);CHKERRQ(ierr);}
1132   PetscFunctionReturn(0);
1133 }
1134 
1135 static PetscErrorCode PCDestroy_BJacobi_Multiproc(PC pc)
1136 {
1137   PC_BJacobi           *jac   = (PC_BJacobi*)pc->data;
1138   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1139   PetscErrorCode       ierr;
1140 
1141   PetscFunctionBegin;
1142   ierr = PCReset_BJacobi_Multiproc(pc);CHKERRQ(ierr);
1143   ierr = KSPDestroy(&jac->ksp[0]);CHKERRQ(ierr);
1144   ierr = PetscFree(jac->ksp);CHKERRQ(ierr);
1145   ierr = PetscSubcommDestroy(&mpjac->psubcomm);CHKERRQ(ierr);
1146 
1147   ierr = PetscFree(mpjac);CHKERRQ(ierr);
1148   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1149   PetscFunctionReturn(0);
1150 }
1151 
1152 static PetscErrorCode PCApply_BJacobi_Multiproc(PC pc,Vec x,Vec y)
1153 {
1154   PC_BJacobi           *jac   = (PC_BJacobi*)pc->data;
1155   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1156   PetscErrorCode       ierr;
1157   PetscScalar          *yarray;
1158   const PetscScalar    *xarray;
1159   KSPConvergedReason   reason;
1160 
1161   PetscFunctionBegin;
1162   /* place x's and y's local arrays into xsub and ysub */
1163   ierr = VecGetArrayRead(x,&xarray);CHKERRQ(ierr);
1164   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
1165   ierr = VecPlaceArray(mpjac->xsub,xarray);CHKERRQ(ierr);
1166   ierr = VecPlaceArray(mpjac->ysub,yarray);CHKERRQ(ierr);
1167 
1168   /* apply preconditioner on each matrix block */
1169   ierr = PetscLogEventBegin(PC_ApplyOnBlocks,jac->ksp[0],mpjac->xsub,mpjac->ysub,0);CHKERRQ(ierr);
1170   ierr = KSPSolve(jac->ksp[0],mpjac->xsub,mpjac->ysub);CHKERRQ(ierr);
1171   ierr = KSPCheckSolve(jac->ksp[0],pc,mpjac->ysub);CHKERRQ(ierr);
1172   ierr = PetscLogEventEnd(PC_ApplyOnBlocks,jac->ksp[0],mpjac->xsub,mpjac->ysub,0);CHKERRQ(ierr);
1173   ierr = KSPGetConvergedReason(jac->ksp[0],&reason);CHKERRQ(ierr);
1174   if (reason == KSP_DIVERGED_PC_FAILED) {
1175     pc->failedreason = PC_SUBPC_ERROR;
1176   }
1177 
1178   ierr = VecResetArray(mpjac->xsub);CHKERRQ(ierr);
1179   ierr = VecResetArray(mpjac->ysub);CHKERRQ(ierr);
1180   ierr = VecRestoreArrayRead(x,&xarray);CHKERRQ(ierr);
1181   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
1182   PetscFunctionReturn(0);
1183 }
1184 
1185 static PetscErrorCode PCMatApply_BJacobi_Multiproc(PC pc,Mat X,Mat Y)
1186 {
1187   PC_BJacobi           *jac   = (PC_BJacobi*)pc->data;
1188   KSPConvergedReason   reason;
1189   Mat                  sX,sY;
1190   const PetscScalar    *x;
1191   PetscScalar          *y;
1192   PetscInt             m,N,lda,ldb;
1193   PetscErrorCode       ierr;
1194 
1195   PetscFunctionBegin;
1196   /* apply preconditioner on each matrix block */
1197   ierr = MatGetLocalSize(X,&m,NULL);CHKERRQ(ierr);
1198   ierr = MatGetSize(X,NULL,&N);CHKERRQ(ierr);
1199   ierr = MatDenseGetLDA(X,&lda);CHKERRQ(ierr);
1200   ierr = MatDenseGetLDA(Y,&ldb);CHKERRQ(ierr);
1201   ierr = MatDenseGetArrayRead(X,&x);CHKERRQ(ierr);
1202   ierr = MatDenseGetArrayWrite(Y,&y);CHKERRQ(ierr);
1203   ierr = MatCreateDense(PetscObjectComm((PetscObject)jac->ksp[0]),m,PETSC_DECIDE,PETSC_DECIDE,N,(PetscScalar*)x,&sX);CHKERRQ(ierr);
1204   ierr = MatCreateDense(PetscObjectComm((PetscObject)jac->ksp[0]),m,PETSC_DECIDE,PETSC_DECIDE,N,y,&sY);CHKERRQ(ierr);
1205   ierr = MatDenseSetLDA(sX,lda);CHKERRQ(ierr);
1206   ierr = MatDenseSetLDA(sY,ldb);CHKERRQ(ierr);
1207   ierr = PetscLogEventBegin(PC_ApplyOnBlocks,jac->ksp[0],X,Y,0);CHKERRQ(ierr);
1208   ierr = KSPMatSolve(jac->ksp[0],sX,sY);CHKERRQ(ierr);
1209   ierr = KSPCheckSolve(jac->ksp[0],pc,NULL);CHKERRQ(ierr);
1210   ierr = PetscLogEventEnd(PC_ApplyOnBlocks,jac->ksp[0],X,Y,0);CHKERRQ(ierr);
1211   ierr = MatDestroy(&sY);CHKERRQ(ierr);
1212   ierr = MatDestroy(&sX);CHKERRQ(ierr);
1213   ierr = MatDenseRestoreArrayWrite(Y,&y);CHKERRQ(ierr);
1214   ierr = MatDenseRestoreArrayRead(X,&x);CHKERRQ(ierr);
1215   ierr = KSPGetConvergedReason(jac->ksp[0],&reason);CHKERRQ(ierr);
1216   if (reason == KSP_DIVERGED_PC_FAILED) {
1217     pc->failedreason = PC_SUBPC_ERROR;
1218   }
1219   PetscFunctionReturn(0);
1220 }
1221 
1222 static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC pc)
1223 {
1224   PC_BJacobi           *jac   = (PC_BJacobi*)pc->data;
1225   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1226   PetscErrorCode       ierr;
1227   PetscInt             m,n;
1228   MPI_Comm             comm,subcomm=0;
1229   const char           *prefix;
1230   PetscBool            wasSetup = PETSC_TRUE;
1231   VecType              vectype;
1232 
1233 
1234   PetscFunctionBegin;
1235   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1236   if (jac->n_local > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only a single block in a subcommunicator is supported");
1237   jac->n_local = 1; /* currently only a single block is supported for a subcommunicator */
1238   if (!pc->setupcalled) {
1239     wasSetup  = PETSC_FALSE;
1240     ierr      = PetscNewLog(pc,&mpjac);CHKERRQ(ierr);
1241     jac->data = (void*)mpjac;
1242 
1243     /* initialize datastructure mpjac */
1244     if (!jac->psubcomm) {
1245       /* Create default contiguous subcommunicatiors if user does not provide them */
1246       ierr = PetscSubcommCreate(comm,&jac->psubcomm);CHKERRQ(ierr);
1247       ierr = PetscSubcommSetNumber(jac->psubcomm,jac->n);CHKERRQ(ierr);
1248       ierr = PetscSubcommSetType(jac->psubcomm,PETSC_SUBCOMM_CONTIGUOUS);CHKERRQ(ierr);
1249       ierr = PetscLogObjectMemory((PetscObject)pc,sizeof(PetscSubcomm));CHKERRQ(ierr);
1250     }
1251     mpjac->psubcomm = jac->psubcomm;
1252     subcomm         = PetscSubcommChild(mpjac->psubcomm);
1253 
1254     /* Get matrix blocks of pmat */
1255     ierr = MatGetMultiProcBlock(pc->pmat,subcomm,MAT_INITIAL_MATRIX,&mpjac->submats);CHKERRQ(ierr);
1256 
1257     /* create a new PC that processors in each subcomm have copy of */
1258     ierr = PetscMalloc1(1,&jac->ksp);CHKERRQ(ierr);
1259     ierr = KSPCreate(subcomm,&jac->ksp[0]);CHKERRQ(ierr);
1260     ierr = KSPSetErrorIfNotConverged(jac->ksp[0],pc->erroriffailure);CHKERRQ(ierr);
1261     ierr = PetscObjectIncrementTabLevel((PetscObject)jac->ksp[0],(PetscObject)pc,1);CHKERRQ(ierr);
1262     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->ksp[0]);CHKERRQ(ierr);
1263     ierr = KSPSetOperators(jac->ksp[0],mpjac->submats,mpjac->submats);CHKERRQ(ierr);
1264     ierr = KSPGetPC(jac->ksp[0],&mpjac->pc);CHKERRQ(ierr);
1265 
1266     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
1267     ierr = KSPSetOptionsPrefix(jac->ksp[0],prefix);CHKERRQ(ierr);
1268     ierr = KSPAppendOptionsPrefix(jac->ksp[0],"sub_");CHKERRQ(ierr);
1269     ierr = KSPGetOptionsPrefix(jac->ksp[0],&prefix);CHKERRQ(ierr);
1270     ierr = MatSetOptionsPrefix(mpjac->submats,prefix);CHKERRQ(ierr);
1271 
1272     /* create dummy vectors xsub and ysub */
1273     ierr = MatGetLocalSize(mpjac->submats,&m,&n);CHKERRQ(ierr);
1274     ierr = VecCreateMPIWithArray(subcomm,1,n,PETSC_DECIDE,NULL,&mpjac->xsub);CHKERRQ(ierr);
1275     ierr = VecCreateMPIWithArray(subcomm,1,m,PETSC_DECIDE,NULL,&mpjac->ysub);CHKERRQ(ierr);
1276     ierr = MatGetVecType(mpjac->submats,&vectype);CHKERRQ(ierr);
1277     ierr = VecSetType(mpjac->xsub,vectype);CHKERRQ(ierr);
1278     ierr = VecSetType(mpjac->ysub,vectype);CHKERRQ(ierr);
1279     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mpjac->xsub);CHKERRQ(ierr);
1280     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mpjac->ysub);CHKERRQ(ierr);
1281 
1282     pc->ops->reset   = PCReset_BJacobi_Multiproc;
1283     pc->ops->destroy = PCDestroy_BJacobi_Multiproc;
1284     pc->ops->apply   = PCApply_BJacobi_Multiproc;
1285     pc->ops->matapply= PCMatApply_BJacobi_Multiproc;
1286   } else { /* pc->setupcalled */
1287     subcomm = PetscSubcommChild(mpjac->psubcomm);
1288     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
1289       /* destroy old matrix blocks, then get new matrix blocks */
1290       if (mpjac->submats) {ierr = MatDestroy(&mpjac->submats);CHKERRQ(ierr);}
1291       ierr = MatGetMultiProcBlock(pc->pmat,subcomm,MAT_INITIAL_MATRIX,&mpjac->submats);CHKERRQ(ierr);
1292     } else {
1293       ierr = MatGetMultiProcBlock(pc->pmat,subcomm,MAT_REUSE_MATRIX,&mpjac->submats);CHKERRQ(ierr);
1294     }
1295     ierr = KSPSetOperators(jac->ksp[0],mpjac->submats,mpjac->submats);CHKERRQ(ierr);
1296   }
1297 
1298   if (!wasSetup && pc->setfromoptionscalled) {
1299     ierr = KSPSetFromOptions(jac->ksp[0]);CHKERRQ(ierr);
1300   }
1301   PetscFunctionReturn(0);
1302 }
1303