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