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