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