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