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