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