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