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