xref: /petsc/src/ksp/pc/impls/bjacobi/bjacobi.c (revision 5a7e1895ba1221176dc63d0bbec6b671c7ff4efe)
1 #define PETSCKSP_DLL
2 
3 /*
4    Defines a block Jacobi preconditioner.
5 */
6 #include "private/pcimpl.h"              /*I "petscpc.h" I*/
7 #include "../src/ksp/pc/impls/bjacobi/bjacobi.h"
8 
9 static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC,Mat,Mat);
10 static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC,Mat,Mat);
11 static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC);
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(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
27   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
28   ierr = MatGetLocalSize(pc->pmat,&M,&N);CHKERRQ(ierr);
29   ierr = MatGetBlockSize(pc->pmat,&bs);CHKERRQ(ierr);
30 
31   if (jac->n > 0 && jac->n < size){
32     ierr = PCSetUp_BJacobi_Multiproc(pc);CHKERRQ(ierr);
33     PetscFunctionReturn(0);
34   }
35 
36   /* --------------------------------------------------------------------------
37       Determines the number of blocks assigned to each processor
38   -----------------------------------------------------------------------------*/
39 
40   /*   local block count  given */
41   if (jac->n_local > 0 && jac->n < 0) {
42     ierr = MPI_Allreduce(&jac->n_local,&jac->n,1,MPIU_INT,MPI_SUM,((PetscObject)pc)->comm);CHKERRQ(ierr);
43     if (jac->l_lens) { /* check that user set these correctly */
44       sum = 0;
45       for (i=0; i<jac->n_local; i++) {
46         if (jac->l_lens[i]/bs*bs !=jac->l_lens[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat blocksize doesn't match block Jacobi layout");
47         sum += jac->l_lens[i];
48       }
49       if (sum != M) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local lens sent incorrectly");
50     } else {
51       ierr = PetscMalloc(jac->n_local*sizeof(PetscInt),&jac->l_lens);CHKERRQ(ierr);
52       for (i=0; i<jac->n_local; i++) {
53         jac->l_lens[i] = bs*((M/bs)/jac->n_local + (((M/bs) % jac->n_local) > i));
54       }
55     }
56   } else if (jac->n > 0 && jac->n_local < 0) { /* global block count given */
57     /* global blocks given: determine which ones are local */
58     if (jac->g_lens) {
59       /* check if the g_lens is has valid entries */
60       for (i=0; i<jac->n; i++) {
61         if (!jac->g_lens[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Zero block not allowed");
62         if (jac->g_lens[i]/bs*bs != jac->g_lens[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat blocksize doesn't match block Jacobi layout");
63       }
64       if (size == 1) {
65         jac->n_local = jac->n;
66         ierr         = PetscMalloc(jac->n_local*sizeof(PetscInt),&jac->l_lens);CHKERRQ(ierr);
67         ierr         = PetscMemcpy(jac->l_lens,jac->g_lens,jac->n_local*sizeof(PetscInt));CHKERRQ(ierr);
68         /* check that user set these correctly */
69         sum = 0;
70         for (i=0; i<jac->n_local; i++) sum += jac->l_lens[i];
71         if (sum != M) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Global lens sent incorrectly");
72       } else {
73         ierr = MatGetOwnershipRange(pc->pmat,&start,&end);CHKERRQ(ierr);
74         /* loop over blocks determing first one owned by me */
75         sum = 0;
76         for (i=0; i<jac->n+1; i++) {
77           if (sum == start) { i_start = i; goto start_1;}
78           if (i < jac->n) sum += jac->g_lens[i];
79         }
80         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Block sizes used in PCBJacobiSetTotalBlocks()\nare 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_COMM_SELF,PETSC_ERR_ARG_SIZ,"Block sizes used in PCBJacobiSetTotalBlocks()\nare not compatible with parallel matrix layout");
87  end_1:
88         jac->n_local = i_end - i_start;
89         ierr         = PetscMalloc(jac->n_local*sizeof(PetscInt),&jac->l_lens);CHKERRQ(ierr);
90         ierr         = PetscMemcpy(jac->l_lens,jac->g_lens+i_start,jac->n_local*sizeof(PetscInt));CHKERRQ(ierr);
91       }
92     } else { /* no global blocks given, determine then using default layout */
93       jac->n_local = jac->n/size + ((jac->n % size) > rank);
94       ierr         = PetscMalloc(jac->n_local*sizeof(PetscInt),&jac->l_lens);CHKERRQ(ierr);
95       for (i=0; i<jac->n_local; i++) {
96         jac->l_lens[i] = ((M/bs)/jac->n_local + (((M/bs) % jac->n_local) > i))*bs;
97         if (!jac->l_lens[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Too many blocks given");
98       }
99     }
100   } else if (jac->n < 0 && jac->n_local < 0) { /* no blocks given */
101     jac->n         = size;
102     jac->n_local   = 1;
103     ierr           = PetscMalloc(sizeof(PetscInt),&jac->l_lens);CHKERRQ(ierr);
104     jac->l_lens[0] = M;
105   }
106   if (jac->n_local < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of blocks is less than number of processors");
107 
108   /* -------------------------
109       Determines mat and pmat
110   ---------------------------*/
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       /* use block from true matrix, not preconditioner matrix for local MatMult() */
121       scall = MAT_INITIAL_MATRIX;
122       if (pc->setupcalled) {
123         if (pc->flag == SAME_NONZERO_PATTERN) {
124           if (jac->tp_mat) {
125             scall = MAT_REUSE_MATRIX;
126             mat   = jac->tp_mat;
127           }
128         } else {
129           if (jac->tp_mat)  {
130             ierr = MatDestroy(jac->tp_mat);CHKERRQ(ierr);
131           }
132         }
133       }
134       if (!f) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"This matrix does not support getting diagonal block");
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         const char *type;
160         ierr = PetscObjectGetType((PetscObject) pc->pmat,&type);CHKERRQ(ierr);
161         SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"This matrix type, %s, does not support getting diagonal block", type);
162       }
163       ierr = (*f)(pc->pmat,&iscopy,scall,&pmat);CHKERRQ(ierr);
164       /* make submatrix have same prefix as entire matrix */
165       ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr);
166       ierr = PetscObjectSetOptionsPrefix((PetscObject)pmat,pprefix);CHKERRQ(ierr);
167       if (iscopy) {
168         jac->tp_pmat = pmat;
169       }
170     } else {
171       pmat = mat;
172     }
173   }
174 
175   /* ------
176      Setup code depends on the number of blocks
177   */
178   if (jac->n_local == 1) {
179     ierr = PCSetUp_BJacobi_Singleblock(pc,mat,pmat);CHKERRQ(ierr);
180   } else {
181     ierr = PCSetUp_BJacobi_Multiblock(pc,mat,pmat);CHKERRQ(ierr);
182   }
183   PetscFunctionReturn(0);
184 }
185 
186 /* Default destroy, if it has never been setup */
187 #undef __FUNCT__
188 #define __FUNCT__ "PCDestroy_BJacobi"
189 static PetscErrorCode PCDestroy_BJacobi(PC pc)
190 {
191   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;
192   PetscErrorCode ierr;
193 
194   PetscFunctionBegin;
195   ierr = PetscFree(jac->g_lens);CHKERRQ(ierr);
196   ierr = PetscFree(jac->l_lens);CHKERRQ(ierr);
197   ierr = PetscFree(jac);CHKERRQ(ierr);
198   PetscFunctionReturn(0);
199 }
200 
201 #undef __FUNCT__
202 #define __FUNCT__ "PCSetFromOptions_BJacobi"
203 
204 static PetscErrorCode PCSetFromOptions_BJacobi(PC pc)
205 {
206   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;
207   PetscErrorCode ierr;
208   PetscInt       blocks;
209   PetscTruth     flg;
210 
211   PetscFunctionBegin;
212   ierr = PetscOptionsHead("Block Jacobi options");CHKERRQ(ierr);
213     ierr = PetscOptionsInt("-pc_bjacobi_blocks","Total number of blocks","PCBJacobiSetTotalBlocks",jac->n,&blocks,&flg);CHKERRQ(ierr);
214     if (flg) {
215       ierr = PCBJacobiSetTotalBlocks(pc,blocks,PETSC_NULL);CHKERRQ(ierr);
216     }
217     flg  = PETSC_FALSE;
218     ierr = PetscOptionsTruth("-pc_bjacobi_truelocal","Use the true matrix, not preconditioner matrix to define matrix vector product in sub-problems","PCBJacobiSetUseTrueLocal",flg,&flg,PETSC_NULL);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,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
239   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&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(((PetscObject)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,((PetscObject)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(((PetscObject)pc)->comm,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(((PetscObject)pc)->comm,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(((PetscObject)pc)->comm,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    Logically 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_CLASSID,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_CLASSID,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 SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG,"Cannot get subsolvers for this preconditioner");
479   PetscFunctionReturn(0);
480 }
481 
482 #undef __FUNCT__
483 #define __FUNCT__ "PCBJacobiSetTotalBlocks"
484 /*@
485    PCBJacobiSetTotalBlocks - Sets the global number of blocks for the block
486    Jacobi preconditioner.
487 
488    Collective on PC
489 
490    Input Parameters:
491 +  pc - the preconditioner context
492 .  blocks - the number of blocks
493 -  lens - [optional] integer array containing the size of each block
494 
495    Options Database Key:
496 .  -pc_bjacobi_blocks <blocks> - Sets the number of global blocks
497 
498    Notes:
499    Currently only a limited number of blocking configurations are supported.
500    All processors sharing the PC must call this routine with the same data.
501 
502    Level: intermediate
503 
504 .keywords:  set, number, Jacobi, global, total, blocks
505 
506 .seealso: PCBJacobiSetUseTrueLocal(), PCBJacobiSetLocalBlocks()
507 @*/
508 PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiSetTotalBlocks(PC pc,PetscInt blocks,const PetscInt lens[])
509 {
510   PetscErrorCode ierr,(*f)(PC,PetscInt,const PetscInt[]);
511 
512   PetscFunctionBegin;
513   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
514   if (blocks <= 0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have positive blocks");
515   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBJacobiSetTotalBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
516   if (f) {
517     ierr = (*f)(pc,blocks,lens);CHKERRQ(ierr);
518   }
519   PetscFunctionReturn(0);
520 }
521 
522 #undef __FUNCT__
523 #define __FUNCT__ "PCBJacobiGetTotalBlocks"
524 /*@C
525    PCBJacobiGetTotalBlocks - Gets the global number of blocks for the block
526    Jacobi preconditioner.
527 
528    Not Collective
529 
530    Input Parameter:
531 .  pc - the preconditioner context
532 
533    Output parameters:
534 +  blocks - the number of blocks
535 -  lens - integer array containing the size of each block
536 
537    Level: intermediate
538 
539 .keywords:  get, number, Jacobi, global, total, blocks
540 
541 .seealso: PCBJacobiSetUseTrueLocal(), PCBJacobiGetLocalBlocks()
542 @*/
543 PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetTotalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
544 {
545   PetscErrorCode ierr,(*f)(PC,PetscInt*, const PetscInt *[]);
546 
547   PetscFunctionBegin;
548   PetscValidHeaderSpecific(pc, PC_CLASSID,1);
549   PetscValidIntPointer(blocks,2);
550   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBJacobiGetTotalBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
551   if (f) {
552     ierr = (*f)(pc,blocks,lens);CHKERRQ(ierr);
553   }
554   PetscFunctionReturn(0);
555 }
556 
557 #undef __FUNCT__
558 #define __FUNCT__ "PCBJacobiSetLocalBlocks"
559 /*@
560    PCBJacobiSetLocalBlocks - Sets the local number of blocks for the block
561    Jacobi preconditioner.
562 
563    Not Collective
564 
565    Input Parameters:
566 +  pc - the preconditioner context
567 .  blocks - the number of blocks
568 -  lens - [optional] integer array containing size of each block
569 
570    Note:
571    Currently only a limited number of blocking configurations are supported.
572 
573    Level: intermediate
574 
575 .keywords: PC, set, number, Jacobi, local, blocks
576 
577 .seealso: PCBJacobiSetUseTrueLocal(), PCBJacobiSetTotalBlocks()
578 @*/
579 PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiSetLocalBlocks(PC pc,PetscInt blocks,const PetscInt lens[])
580 {
581   PetscErrorCode ierr,(*f)(PC,PetscInt,const PetscInt []);
582 
583   PetscFunctionBegin;
584   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
585   if (blocks < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have nonegative blocks");
586   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBJacobiSetLocalBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
587   if (f) {
588     ierr = (*f)(pc,blocks,lens);CHKERRQ(ierr);
589   }
590   PetscFunctionReturn(0);
591 }
592 
593 #undef __FUNCT__
594 #define __FUNCT__ "PCBJacobiGetLocalBlocks"
595 /*@C
596    PCBJacobiGetLocalBlocks - Gets the local number of blocks for the block
597    Jacobi preconditioner.
598 
599    Not Collective
600 
601    Input Parameters:
602 +  pc - the preconditioner context
603 .  blocks - the number of blocks
604 -  lens - [optional] integer array containing size of each block
605 
606    Note:
607    Currently only a limited number of blocking configurations are supported.
608 
609    Level: intermediate
610 
611 .keywords: PC, get, number, Jacobi, local, blocks
612 
613 .seealso: PCBJacobiSetUseTrueLocal(), PCBJacobiGetTotalBlocks()
614 @*/
615 PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetLocalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
616 {
617   PetscErrorCode ierr,(*f)(PC,PetscInt*, const PetscInt *[]);
618 
619   PetscFunctionBegin;
620   PetscValidHeaderSpecific(pc, PC_CLASSID,1);
621   PetscValidIntPointer(blocks,2);
622   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBJacobiGetLocalBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
623   if (f) {
624     ierr = (*f)(pc,blocks,lens);CHKERRQ(ierr);
625   }
626   PetscFunctionReturn(0);
627 }
628 
629 /* -----------------------------------------------------------------------------------*/
630 
631 /*MC
632    PCBJACOBI - Use block Jacobi preconditioning, each block is (approximately) solved with
633            its own KSP object.
634 
635    Options Database Keys:
636 .  -pc_bjacobi_truelocal - Activates PCBJacobiSetUseTrueLocal()
637 
638    Notes: Each processor can have one or more blocks, but a block cannot be shared by more
639      than one processor. Defaults to one block per processor.
640 
641      To set options on the solvers for each block append -sub_ to all the KSP, KSP, and PC
642         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
643 
644      To set the options on the solvers separate for each block call PCBJacobiGetSubKSP()
645          and set the options directly on the resulting KSP object (you can access its PC
646          KSPGetPC())
647 
648    Level: beginner
649 
650    Concepts: block Jacobi
651 
652 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
653            PCASM, PCBJacobiSetUseTrueLocal(), PCBJacobiGetSubKSP(), PCBJacobiSetTotalBlocks(),
654            PCBJacobiSetLocalBlocks(), PCSetModifySubmatrices()
655 M*/
656 
657 EXTERN_C_BEGIN
658 #undef __FUNCT__
659 #define __FUNCT__ "PCCreate_BJacobi"
660 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_BJacobi(PC pc)
661 {
662   PetscErrorCode ierr;
663   PetscMPIInt    rank;
664   PC_BJacobi     *jac;
665 
666   PetscFunctionBegin;
667   ierr = PetscNewLog(pc,PC_BJacobi,&jac);CHKERRQ(ierr);
668   ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
669   pc->ops->apply              = 0;
670   pc->ops->applytranspose     = 0;
671   pc->ops->setup              = PCSetUp_BJacobi;
672   pc->ops->destroy            = PCDestroy_BJacobi;
673   pc->ops->setfromoptions     = PCSetFromOptions_BJacobi;
674   pc->ops->view               = PCView_BJacobi;
675   pc->ops->applyrichardson    = 0;
676 
677   pc->data               = (void*)jac;
678   jac->n                 = -1;
679   jac->n_local           = -1;
680   jac->first_local       = rank;
681   jac->ksp              = 0;
682   jac->use_true_local    = PETSC_FALSE;
683   jac->same_local_solves = PETSC_TRUE;
684   jac->g_lens            = 0;
685   jac->l_lens            = 0;
686   jac->tp_mat            = 0;
687   jac->tp_pmat           = 0;
688 
689   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiSetUseTrueLocal_C",
690                     "PCBJacobiSetUseTrueLocal_BJacobi",
691                     PCBJacobiSetUseTrueLocal_BJacobi);CHKERRQ(ierr);
692   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiGetSubKSP_C","PCBJacobiGetSubKSP_BJacobi",
693                     PCBJacobiGetSubKSP_BJacobi);CHKERRQ(ierr);
694   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiSetTotalBlocks_C","PCBJacobiSetTotalBlocks_BJacobi",
695                     PCBJacobiSetTotalBlocks_BJacobi);CHKERRQ(ierr);
696   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiGetTotalBlocks_C","PCBJacobiGetTotalBlocks_BJacobi",
697                     PCBJacobiGetTotalBlocks_BJacobi);CHKERRQ(ierr);
698   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiSetLocalBlocks_C","PCBJacobiSetLocalBlocks_BJacobi",
699                     PCBJacobiSetLocalBlocks_BJacobi);CHKERRQ(ierr);
700   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiGetLocalBlocks_C","PCBJacobiGetLocalBlocks_BJacobi",
701                     PCBJacobiGetLocalBlocks_BJacobi);CHKERRQ(ierr);
702 
703   PetscFunctionReturn(0);
704 }
705 EXTERN_C_END
706 
707 /* --------------------------------------------------------------------------------------------*/
708 /*
709         These are for a single block per processor; works for AIJ, BAIJ; Seq and MPI
710 */
711 #undef __FUNCT__
712 #define __FUNCT__ "PCDestroy_BJacobi_Singleblock"
713 PetscErrorCode PCDestroy_BJacobi_Singleblock(PC pc)
714 {
715   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
716   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
717   PetscErrorCode         ierr;
718 
719   PetscFunctionBegin;
720   /*
721         If the on processor block had to be generated via a MatGetDiagonalBlock()
722      that creates a copy, this frees the space
723   */
724   if (jac->tp_mat) {
725     ierr = MatDestroy(jac->tp_mat);CHKERRQ(ierr);
726   }
727   if (jac->tp_pmat) {
728     ierr = MatDestroy(jac->tp_pmat);CHKERRQ(ierr);
729   }
730 
731   ierr = KSPDestroy(jac->ksp[0]);CHKERRQ(ierr);
732   ierr = PetscFree(jac->ksp);CHKERRQ(ierr);
733   ierr = VecDestroy(bjac->x);CHKERRQ(ierr);
734   ierr = VecDestroy(bjac->y);CHKERRQ(ierr);
735   ierr = PetscFree(jac->l_lens);CHKERRQ(ierr);
736   ierr = PetscFree(jac->g_lens);CHKERRQ(ierr);
737   ierr = PetscFree(bjac);CHKERRQ(ierr);
738   ierr = PetscFree(jac);CHKERRQ(ierr);
739   PetscFunctionReturn(0);
740 }
741 
742 #undef __FUNCT__
743 #define __FUNCT__ "PCSetUpOnBlocks_BJacobi_Singleblock"
744 PetscErrorCode PCSetUpOnBlocks_BJacobi_Singleblock(PC pc)
745 {
746   PetscErrorCode ierr;
747   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;
748 
749   PetscFunctionBegin;
750   ierr = KSPSetUp(jac->ksp[0]);CHKERRQ(ierr);
751   PetscFunctionReturn(0);
752 }
753 
754 #undef __FUNCT__
755 #define __FUNCT__ "PCApply_BJacobi_Singleblock"
756 PetscErrorCode PCApply_BJacobi_Singleblock(PC pc,Vec x,Vec y)
757 {
758   PetscErrorCode         ierr;
759   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
760   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
761   PetscScalar            *x_array,*y_array;
762 
763   PetscFunctionBegin;
764   /*
765       The VecPlaceArray() is to avoid having to copy the
766     y vector into the bjac->x vector. The reason for
767     the bjac->x vector is that we need a sequential vector
768     for the sequential solve.
769   */
770   ierr = VecGetArray(x,&x_array);CHKERRQ(ierr);
771   ierr = VecGetArray(y,&y_array);CHKERRQ(ierr);
772   ierr = VecPlaceArray(bjac->x,x_array);CHKERRQ(ierr);
773   ierr = VecPlaceArray(bjac->y,y_array);CHKERRQ(ierr);
774   ierr = KSPSolve(jac->ksp[0],bjac->x,bjac->y);CHKERRQ(ierr);
775   ierr = VecResetArray(bjac->x);CHKERRQ(ierr);
776   ierr = VecResetArray(bjac->y);CHKERRQ(ierr);
777   ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr);
778   ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr);
779   PetscFunctionReturn(0);
780 }
781 
782 #undef __FUNCT__
783 #define __FUNCT__ "PCApplySymmetricLeft_BJacobi_Singleblock"
784 PetscErrorCode PCApplySymmetricLeft_BJacobi_Singleblock(PC pc,Vec x,Vec y)
785 {
786   PetscErrorCode         ierr;
787   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
788   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
789   PetscScalar            *x_array,*y_array;
790   PC                     subpc;
791 
792   PetscFunctionBegin;
793   /*
794       The VecPlaceArray() is to avoid having to copy the
795     y vector into the bjac->x vector. The reason for
796     the bjac->x vector is that we need a sequential vector
797     for the sequential solve.
798   */
799   ierr = VecGetArray(x,&x_array);CHKERRQ(ierr);
800   ierr = VecGetArray(y,&y_array);CHKERRQ(ierr);
801   ierr = VecPlaceArray(bjac->x,x_array);CHKERRQ(ierr);
802   ierr = VecPlaceArray(bjac->y,y_array);CHKERRQ(ierr);
803 
804   /* apply the symmetric left portion of the inner PC operator */
805   /* note this by-passes the inner KSP and its options completely */
806 
807   ierr = KSPGetPC(jac->ksp[0],&subpc);CHKERRQ(ierr);
808   ierr = PCApplySymmetricLeft(subpc,bjac->x,bjac->y);CHKERRQ(ierr);
809   ierr = VecResetArray(bjac->x);CHKERRQ(ierr);
810   ierr = VecResetArray(bjac->y);CHKERRQ(ierr);
811 
812   ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr);
813   ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr);
814   PetscFunctionReturn(0);
815 }
816 
817 #undef __FUNCT__
818 #define __FUNCT__ "PCApplySymmetricRight_BJacobi_Singleblock"
819 PetscErrorCode PCApplySymmetricRight_BJacobi_Singleblock(PC pc,Vec x,Vec y)
820 {
821   PetscErrorCode         ierr;
822   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
823   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
824   PetscScalar            *x_array,*y_array;
825   PC                     subpc;
826 
827   PetscFunctionBegin;
828   /*
829       The VecPlaceArray() is to avoid having to copy the
830     y vector into the bjac->x vector. The reason for
831     the bjac->x vector is that we need a sequential vector
832     for the sequential solve.
833   */
834   ierr = VecGetArray(x,&x_array);CHKERRQ(ierr);
835   ierr = VecGetArray(y,&y_array);CHKERRQ(ierr);
836   ierr = VecPlaceArray(bjac->x,x_array);CHKERRQ(ierr);
837   ierr = VecPlaceArray(bjac->y,y_array);CHKERRQ(ierr);
838 
839   /* apply the symmetric right portion of the inner PC operator */
840   /* note this by-passes the inner KSP and its options completely */
841 
842   ierr = KSPGetPC(jac->ksp[0],&subpc);CHKERRQ(ierr);
843   ierr = PCApplySymmetricRight(subpc,bjac->x,bjac->y);CHKERRQ(ierr);
844 
845   ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr);
846   ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr);
847   PetscFunctionReturn(0);
848 }
849 
850 #undef __FUNCT__
851 #define __FUNCT__ "PCApplyTranspose_BJacobi_Singleblock"
852 PetscErrorCode PCApplyTranspose_BJacobi_Singleblock(PC pc,Vec x,Vec y)
853 {
854   PetscErrorCode         ierr;
855   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
856   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
857   PetscScalar            *x_array,*y_array;
858 
859   PetscFunctionBegin;
860   /*
861       The VecPlaceArray() is to avoid having to copy the
862     y vector into the bjac->x vector. The reason for
863     the bjac->x vector is that we need a sequential vector
864     for the sequential solve.
865   */
866   ierr = VecGetArray(x,&x_array);CHKERRQ(ierr);
867   ierr = VecGetArray(y,&y_array);CHKERRQ(ierr);
868   ierr = VecPlaceArray(bjac->x,x_array);CHKERRQ(ierr);
869   ierr = VecPlaceArray(bjac->y,y_array);CHKERRQ(ierr);
870   ierr = KSPSolveTranspose(jac->ksp[0],bjac->x,bjac->y);CHKERRQ(ierr);
871   ierr = VecResetArray(bjac->x);CHKERRQ(ierr);
872   ierr = VecResetArray(bjac->y);CHKERRQ(ierr);
873   ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr);
874   ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr);
875   PetscFunctionReturn(0);
876 }
877 
878 #undef __FUNCT__
879 #define __FUNCT__ "PCSetUp_BJacobi_Singleblock"
880 static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC pc,Mat mat,Mat pmat)
881 {
882   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
883   PetscErrorCode         ierr;
884   PetscInt               m;
885   KSP                    ksp;
886   Vec                    x,y;
887   PC_BJacobi_Singleblock *bjac;
888   PetscTruth             wasSetup;
889 
890   PetscFunctionBegin;
891 
892   /* set default direct solver with no Krylov method */
893   if (!pc->setupcalled) {
894     const char *prefix;
895     wasSetup = PETSC_FALSE;
896     ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr);
897     ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
898     ierr = PetscLogObjectParent(pc,ksp);CHKERRQ(ierr);
899     ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
900     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
901     ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
902     ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
903     /*
904       The reason we need to generate these vectors is to serve
905       as the right-hand side and solution vector for the solve on the
906       block. We do not need to allocate space for the vectors since
907       that is provided via VecPlaceArray() just before the call to
908       KSPSolve() on the block.
909     */
910     ierr = MatGetSize(pmat,&m,&m);CHKERRQ(ierr);
911     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,m,PETSC_NULL,&x);CHKERRQ(ierr);
912     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,m,PETSC_NULL,&y);CHKERRQ(ierr);
913     ierr = PetscLogObjectParent(pc,x);CHKERRQ(ierr);
914     ierr = PetscLogObjectParent(pc,y);CHKERRQ(ierr);
915 
916     pc->ops->destroy             = PCDestroy_BJacobi_Singleblock;
917     pc->ops->apply               = PCApply_BJacobi_Singleblock;
918     pc->ops->applysymmetricleft  = PCApplySymmetricLeft_BJacobi_Singleblock;
919     pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Singleblock;
920     pc->ops->applytranspose      = PCApplyTranspose_BJacobi_Singleblock;
921     pc->ops->setuponblocks       = PCSetUpOnBlocks_BJacobi_Singleblock;
922 
923     ierr = PetscNewLog(pc,PC_BJacobi_Singleblock,&bjac);CHKERRQ(ierr);
924     bjac->x      = x;
925     bjac->y      = y;
926 
927     ierr = PetscMalloc(sizeof(KSP),&jac->ksp);CHKERRQ(ierr);
928     jac->ksp[0] = ksp;
929     jac->data    = (void*)bjac;
930   } else {
931     wasSetup = PETSC_TRUE;
932     ksp = jac->ksp[0];
933     bjac = (PC_BJacobi_Singleblock *)jac->data;
934   }
935   if (jac->use_true_local) {
936     ierr = KSPSetOperators(ksp,mat,pmat,pc->flag);CHKERRQ(ierr);
937   }  else {
938     ierr = KSPSetOperators(ksp,pmat,pmat,pc->flag);CHKERRQ(ierr);
939   }
940   if (!wasSetup && pc->setfromoptionscalled) {
941     ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
942   }
943   PetscFunctionReturn(0);
944 }
945 
946 /* ---------------------------------------------------------------------------------------------*/
947 
948 #undef __FUNCT__
949 #define __FUNCT__ "PCDestroy_BJacobi_Multiblock"
950 PetscErrorCode PCDestroy_BJacobi_Multiblock(PC pc)
951 {
952   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
953   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
954   PetscErrorCode        ierr;
955   PetscInt              i;
956 
957   PetscFunctionBegin;
958   ierr = MatDestroyMatrices(jac->n_local,&bjac->pmat);CHKERRQ(ierr);
959   if (jac->use_true_local) {
960     ierr = MatDestroyMatrices(jac->n_local,&bjac->mat);CHKERRQ(ierr);
961   }
962 
963   /*
964         If the on processor block had to be generated via a MatGetDiagonalBlock()
965      that creates a copy, this frees the space
966   */
967   if (jac->tp_mat) {
968     ierr = MatDestroy(jac->tp_mat);CHKERRQ(ierr);
969   }
970   if (jac->tp_pmat) {
971     ierr = MatDestroy(jac->tp_pmat);CHKERRQ(ierr);
972   }
973 
974   for (i=0; i<jac->n_local; i++) {
975     ierr = KSPDestroy(jac->ksp[i]);CHKERRQ(ierr);
976     ierr = VecDestroy(bjac->x[i]);CHKERRQ(ierr);
977     ierr = VecDestroy(bjac->y[i]);CHKERRQ(ierr);
978     ierr = ISDestroy(bjac->is[i]);CHKERRQ(ierr);
979   }
980   ierr = PetscFree(jac->ksp);CHKERRQ(ierr);
981   ierr = PetscFree2(bjac->x,bjac->y);CHKERRQ(ierr);
982   ierr = PetscFree(bjac->starts);CHKERRQ(ierr);
983   ierr = PetscFree(bjac->is);CHKERRQ(ierr);
984   ierr = PetscFree(bjac);CHKERRQ(ierr);
985   ierr = PetscFree(jac->l_lens);CHKERRQ(ierr);
986   ierr = PetscFree(jac->g_lens);CHKERRQ(ierr);
987   ierr = PetscFree(jac);CHKERRQ(ierr);
988   PetscFunctionReturn(0);
989 }
990 
991 #undef __FUNCT__
992 #define __FUNCT__ "PCSetUpOnBlocks_BJacobi_Multiblock"
993 PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiblock(PC pc)
994 {
995   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;
996   PetscErrorCode ierr;
997   PetscInt       i,n_local = jac->n_local;
998 
999   PetscFunctionBegin;
1000   for (i=0; i<n_local; i++) {
1001     ierr = KSPSetUp(jac->ksp[i]);CHKERRQ(ierr);
1002   }
1003   PetscFunctionReturn(0);
1004 }
1005 
1006 /*
1007       Preconditioner for block Jacobi
1008 */
1009 #undef __FUNCT__
1010 #define __FUNCT__ "PCApply_BJacobi_Multiblock"
1011 PetscErrorCode PCApply_BJacobi_Multiblock(PC pc,Vec x,Vec y)
1012 {
1013   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
1014   PetscErrorCode        ierr;
1015   PetscInt              i,n_local = jac->n_local;
1016   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
1017   PetscScalar           *xin,*yin;
1018 
1019   PetscFunctionBegin;
1020   ierr = VecGetArray(x,&xin);CHKERRQ(ierr);
1021   ierr = VecGetArray(y,&yin);CHKERRQ(ierr);
1022   for (i=0; i<n_local; i++) {
1023     /*
1024        To avoid copying the subvector from x into a workspace we instead
1025        make the workspace vector array point to the subpart of the array of
1026        the global vector.
1027     */
1028     ierr = VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);CHKERRQ(ierr);
1029     ierr = VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);CHKERRQ(ierr);
1030 
1031     ierr = PetscLogEventBegin(PC_SetUpOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);CHKERRQ(ierr);
1032     ierr = KSPSolve(jac->ksp[i],bjac->x[i],bjac->y[i]);CHKERRQ(ierr);
1033     ierr = PetscLogEventEnd(PC_SetUpOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);CHKERRQ(ierr);
1034 
1035     ierr = VecResetArray(bjac->x[i]);CHKERRQ(ierr);
1036     ierr = VecResetArray(bjac->y[i]);CHKERRQ(ierr);
1037   }
1038   ierr = VecRestoreArray(x,&xin);CHKERRQ(ierr);
1039   ierr = VecRestoreArray(y,&yin);CHKERRQ(ierr);
1040   PetscFunctionReturn(0);
1041 }
1042 
1043 /*
1044       Preconditioner for block Jacobi
1045 */
1046 #undef __FUNCT__
1047 #define __FUNCT__ "PCApplyTranspose_BJacobi_Multiblock"
1048 PetscErrorCode PCApplyTranspose_BJacobi_Multiblock(PC pc,Vec x,Vec y)
1049 {
1050   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
1051   PetscErrorCode        ierr;
1052   PetscInt              i,n_local = jac->n_local;
1053   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
1054   PetscScalar           *xin,*yin;
1055 
1056   PetscFunctionBegin;
1057   ierr = VecGetArray(x,&xin);CHKERRQ(ierr);
1058   ierr = VecGetArray(y,&yin);CHKERRQ(ierr);
1059   for (i=0; i<n_local; i++) {
1060     /*
1061        To avoid copying the subvector from x into a workspace we instead
1062        make the workspace vector array point to the subpart of the array of
1063        the global vector.
1064     */
1065     ierr = VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);CHKERRQ(ierr);
1066     ierr = VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);CHKERRQ(ierr);
1067 
1068     ierr = PetscLogEventBegin(PC_ApplyTransposeOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);CHKERRQ(ierr);
1069     ierr = KSPSolveTranspose(jac->ksp[i],bjac->x[i],bjac->y[i]);CHKERRQ(ierr);
1070     ierr = PetscLogEventEnd(PC_ApplyTransposeOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);CHKERRQ(ierr);
1071   }
1072   ierr = VecRestoreArray(x,&xin);CHKERRQ(ierr);
1073   ierr = VecRestoreArray(y,&yin);CHKERRQ(ierr);
1074   PetscFunctionReturn(0);
1075 }
1076 
1077 #undef __FUNCT__
1078 #define __FUNCT__ "PCSetUp_BJacobi_Multiblock"
1079 static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC pc,Mat mat,Mat pmat)
1080 {
1081   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
1082   PetscErrorCode         ierr;
1083   PetscInt               m,n_local,N,M,start,i;
1084   const char             *prefix,*pprefix,*mprefix;
1085   KSP                    ksp;
1086   Vec                    x,y;
1087   PC_BJacobi_Multiblock  *bjac = (PC_BJacobi_Multiblock*)jac->data;
1088   PC                     subpc;
1089   IS                     is;
1090   MatReuse               scall = MAT_REUSE_MATRIX;
1091 
1092   PetscFunctionBegin;
1093   ierr = MatGetLocalSize(pc->pmat,&M,&N);CHKERRQ(ierr);
1094 
1095   n_local = jac->n_local;
1096 
1097   if (jac->use_true_local) {
1098     PetscTruth same;
1099     ierr = PetscTypeCompare((PetscObject)mat,((PetscObject)pmat)->type_name,&same);CHKERRQ(ierr);
1100     if (!same) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Matrices not of same type");
1101   }
1102 
1103   if (!pc->setupcalled) {
1104     scall                  = MAT_INITIAL_MATRIX;
1105     pc->ops->destroy       = PCDestroy_BJacobi_Multiblock;
1106     pc->ops->apply         = PCApply_BJacobi_Multiblock;
1107     pc->ops->applytranspose= PCApplyTranspose_BJacobi_Multiblock;
1108     pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiblock;
1109 
1110     ierr = PetscNewLog(pc,PC_BJacobi_Multiblock,&bjac);CHKERRQ(ierr);
1111     ierr = PetscMalloc(n_local*sizeof(KSP),&jac->ksp);CHKERRQ(ierr);
1112     ierr = PetscLogObjectMemory(pc,sizeof(n_local*sizeof(KSP)));CHKERRQ(ierr);
1113     ierr = PetscMalloc2(n_local,Vec,&bjac->x,n_local,Vec,&bjac->y);CHKERRQ(ierr);
1114     ierr = PetscMalloc(n_local*sizeof(PetscScalar),&bjac->starts);CHKERRQ(ierr);
1115     ierr = PetscLogObjectMemory(pc,sizeof(n_local*sizeof(PetscScalar)));CHKERRQ(ierr);
1116 
1117     jac->data    = (void*)bjac;
1118     ierr = PetscMalloc(n_local*sizeof(IS),&bjac->is);CHKERRQ(ierr);
1119     ierr = PetscLogObjectMemory(pc,sizeof(n_local*sizeof(IS)));CHKERRQ(ierr);
1120 
1121     start = 0;
1122     for (i=0; i<n_local; i++) {
1123       ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr);
1124       ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1125       ierr = PetscLogObjectParent(pc,ksp);CHKERRQ(ierr);
1126       ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
1127       ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr);
1128       ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
1129       ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
1130       ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
1131 
1132       m = jac->l_lens[i];
1133 
1134       /*
1135       The reason we need to generate these vectors is to serve
1136       as the right-hand side and solution vector for the solve on the
1137       block. We do not need to allocate space for the vectors since
1138       that is provided via VecPlaceArray() just before the call to
1139       KSPSolve() on the block.
1140 
1141       */
1142       ierr = VecCreateSeq(PETSC_COMM_SELF,m,&x);CHKERRQ(ierr);
1143       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,m,PETSC_NULL,&y);CHKERRQ(ierr);
1144       ierr = PetscLogObjectParent(pc,x);CHKERRQ(ierr);
1145       ierr = PetscLogObjectParent(pc,y);CHKERRQ(ierr);
1146       bjac->x[i]      = x;
1147       bjac->y[i]      = y;
1148       bjac->starts[i] = start;
1149       jac->ksp[i]    = ksp;
1150 
1151       ierr = ISCreateStride(PETSC_COMM_SELF,m,start,1,&is);CHKERRQ(ierr);
1152       bjac->is[i] = is;
1153       ierr = PetscLogObjectParent(pc,is);CHKERRQ(ierr);
1154 
1155       start += m;
1156     }
1157   } else {
1158     bjac = (PC_BJacobi_Multiblock*)jac->data;
1159     /*
1160        Destroy the blocks from the previous iteration
1161     */
1162     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
1163       ierr = MatDestroyMatrices(n_local,&bjac->pmat);CHKERRQ(ierr);
1164       if (jac->use_true_local) {
1165         ierr = MatDestroyMatrices(n_local,&bjac->mat);CHKERRQ(ierr);
1166       }
1167       scall = MAT_INITIAL_MATRIX;
1168     }
1169   }
1170 
1171   ierr = MatGetSubMatrices(pmat,n_local,bjac->is,bjac->is,scall,&bjac->pmat);CHKERRQ(ierr);
1172   if (jac->use_true_local) {
1173     ierr = PetscObjectGetOptionsPrefix((PetscObject)mat,&mprefix);CHKERRQ(ierr);
1174     ierr = MatGetSubMatrices(mat,n_local,bjac->is,bjac->is,scall,&bjac->mat);CHKERRQ(ierr);
1175   }
1176   /* Return control to the user so that the submatrices can be modified (e.g., to apply
1177      different boundary conditions for the submatrices than for the global problem) */
1178   ierr = PCModifySubMatrices(pc,n_local,bjac->is,bjac->is,bjac->pmat,pc->modifysubmatricesP);CHKERRQ(ierr);
1179 
1180   ierr = PetscObjectGetOptionsPrefix((PetscObject)pmat,&pprefix);CHKERRQ(ierr);
1181   for (i=0; i<n_local; i++) {
1182     ierr = PetscLogObjectParent(pc,bjac->pmat[i]);CHKERRQ(ierr);
1183     ierr = PetscObjectSetOptionsPrefix((PetscObject)bjac->pmat[i],pprefix);CHKERRQ(ierr);
1184     if (jac->use_true_local) {
1185       ierr = PetscLogObjectParent(pc,bjac->mat[i]);CHKERRQ(ierr);
1186       ierr = PetscObjectSetOptionsPrefix((PetscObject)bjac->mat[i],mprefix);CHKERRQ(ierr);
1187       ierr = KSPSetOperators(jac->ksp[i],bjac->mat[i],bjac->pmat[i],pc->flag);CHKERRQ(ierr);
1188     } else {
1189       ierr = KSPSetOperators(jac->ksp[i],bjac->pmat[i],bjac->pmat[i],pc->flag);CHKERRQ(ierr);
1190     }
1191     if(pc->setfromoptionscalled){
1192       ierr = KSPSetFromOptions(jac->ksp[i]);CHKERRQ(ierr);
1193     }
1194   }
1195   PetscFunctionReturn(0);
1196 }
1197 
1198 /* ---------------------------------------------------------------------------------------------*/
1199 /*
1200       These are for a single block with multiple processes;
1201 */
1202 #undef __FUNCT__
1203 #define __FUNCT__ "PCDestroy_BJacobi_Multiproc"
1204 static PetscErrorCode PCDestroy_BJacobi_Multiproc(PC pc)
1205 {
1206   PC_BJacobi           *jac = (PC_BJacobi*)pc->data;
1207   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1208   PetscErrorCode       ierr;
1209 
1210   PetscFunctionBegin;
1211   if (mpjac->ysub){ierr = VecDestroy(mpjac->ysub);CHKERRQ(ierr);}
1212   if (mpjac->xsub){ierr = VecDestroy(mpjac->xsub);CHKERRQ(ierr);}
1213   if (mpjac->submats){ierr = MatDestroy(mpjac->submats);CHKERRQ(ierr);}
1214   if (mpjac->ksp){ierr = KSPDestroy(mpjac->ksp);CHKERRQ(ierr);}
1215   if (mpjac->psubcomm){ierr = PetscSubcommDestroy(mpjac->psubcomm);CHKERRQ(ierr);}
1216 
1217   ierr = PetscFree(mpjac);CHKERRQ(ierr);
1218   ierr = PetscFree(jac);CHKERRQ(ierr);
1219   PetscFunctionReturn(0);
1220 }
1221 
1222 #undef __FUNCT__
1223 #define __FUNCT__ "PCApply_BJacobi_Multiproc"
1224 static PetscErrorCode PCApply_BJacobi_Multiproc(PC pc,Vec x,Vec y)
1225 {
1226   PC_BJacobi           *jac = (PC_BJacobi*)pc->data;
1227   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1228   PetscErrorCode       ierr;
1229   PetscScalar          *xarray,*yarray;
1230 
1231   PetscFunctionBegin;
1232   /* place x's and y's local arrays into xsub and ysub */
1233   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
1234   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
1235   ierr = VecPlaceArray(mpjac->xsub,xarray);CHKERRQ(ierr);
1236   ierr = VecPlaceArray(mpjac->ysub,yarray);CHKERRQ(ierr);
1237 
1238   /* apply preconditioner on each matrix block */
1239   ierr = KSPSolve(mpjac->ksp,mpjac->xsub,mpjac->ysub);CHKERRQ(ierr);
1240 
1241   ierr = VecResetArray(mpjac->xsub);CHKERRQ(ierr);
1242   ierr = VecResetArray(mpjac->ysub);CHKERRQ(ierr);
1243   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
1244   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
1245   PetscFunctionReturn(0);
1246 }
1247 
1248 extern PetscErrorCode MatGetMultiProcBlock_MPIAIJ(Mat,MPI_Comm,Mat*);
1249 #undef __FUNCT__
1250 #define __FUNCT__ "PCSetUp_BJacobi_Multiproc"
1251 static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC pc)
1252 {
1253   PC_BJacobi           *jac = (PC_BJacobi*)pc->data;
1254   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1255   PetscErrorCode       ierr;
1256   PetscInt             m,n;
1257   MPI_Comm             comm = ((PetscObject)pc)->comm,subcomm;
1258   const char           *prefix;
1259 
1260   PetscFunctionBegin;
1261   if (jac->n_local > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only a single block in a subcommunicator is supported");
1262   jac->n_local = 1; /* currently only a single block is supported for a subcommunicator */
1263   if (!pc->setupcalled) {
1264     ierr = PetscNewLog(pc,PC_BJacobi_Multiproc,&mpjac);CHKERRQ(ierr);
1265     jac->data = (void*)mpjac;
1266 
1267     /* initialize datastructure mpjac */
1268     mpjac->nsubcomm = jac->n;
1269 
1270     if (!jac->psubcomm) {
1271       /* Create default contiguous subcommunicatiors if user does not provide them */
1272       ierr = PetscSubcommCreate(comm,mpjac->nsubcomm,PETSC_SUBCOMM_CONTIGUOUS,&jac->psubcomm);CHKERRQ(ierr);
1273       ierr = PetscLogObjectMemory(pc,sizeof(PetscSubcomm));CHKERRQ(ierr);
1274     }
1275     mpjac->psubcomm = jac->psubcomm;
1276     subcomm         = mpjac->psubcomm->comm;
1277 
1278     /* Get matrix blocks of pmat */
1279     ierr = MatGetMultiProcBlock_MPIAIJ(pc->pmat,subcomm,&mpjac->submats);CHKERRQ(ierr);
1280 
1281     /* create a new PC that processors in each subcomm have copy of */
1282     ierr = KSPCreate(subcomm,&mpjac->ksp);CHKERRQ(ierr);
1283     ierr = PetscObjectIncrementTabLevel((PetscObject)mpjac->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1284     ierr = PetscLogObjectParent(pc,mpjac->ksp);CHKERRQ(ierr);
1285     ierr = KSPSetOperators(mpjac->ksp,mpjac->submats,mpjac->submats,pc->flag);CHKERRQ(ierr);
1286     ierr = KSPGetPC(mpjac->ksp,&mpjac->pc);CHKERRQ(ierr);
1287 
1288     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
1289     ierr = KSPSetOptionsPrefix(mpjac->ksp,prefix);CHKERRQ(ierr);
1290     ierr = KSPAppendOptionsPrefix(mpjac->ksp,"sub_");CHKERRQ(ierr);
1291     /*
1292       PetscMPIInt rank,subsize,subrank;
1293       ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1294       ierr = MPI_Comm_size(subcomm,&subsize);CHKERRQ(ierr);
1295       ierr = MPI_Comm_rank(subcomm,&subrank);CHKERRQ(ierr);
1296 
1297       ierr = MatGetLocalSize(mpjac->submats,&m,PETSC_NULL);CHKERRQ(ierr);
1298       ierr = MatGetSize(mpjac->submats,&n,PETSC_NULL);CHKERRQ(ierr);
1299       ierr = PetscSynchronizedPrintf(comm,"[%d], sub-size %d,sub-rank %d\n",rank,subsize,subrank);
1300       ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr);
1301     */
1302 
1303     /* create dummy vectors xsub and ysub */
1304     ierr = MatGetLocalSize(mpjac->submats,&m,&n);CHKERRQ(ierr);
1305     ierr = VecCreateMPIWithArray(subcomm,n,PETSC_DECIDE,PETSC_NULL,&mpjac->xsub);CHKERRQ(ierr);
1306     ierr = VecCreateMPIWithArray(subcomm,m,PETSC_DECIDE,PETSC_NULL,&mpjac->ysub);CHKERRQ(ierr);
1307     ierr = PetscLogObjectParent(pc,mpjac->xsub);CHKERRQ(ierr);
1308     ierr = PetscLogObjectParent(pc,mpjac->ysub);CHKERRQ(ierr);
1309 
1310     pc->ops->destroy = PCDestroy_BJacobi_Multiproc;
1311     pc->ops->apply   = PCApply_BJacobi_Multiproc;
1312   }
1313 
1314   if (pc->setupcalled && pc->flag == DIFFERENT_NONZERO_PATTERN) {
1315     /* destroy old matrix blocks, then get new matrix blocks */
1316     if (mpjac->submats) {
1317       ierr = MatDestroy(mpjac->submats);CHKERRQ(ierr);
1318       ierr = MatGetMultiProcBlock_MPIAIJ(pc->pmat,subcomm,&mpjac->submats);CHKERRQ(ierr);
1319       ierr = KSPSetOperators(mpjac->ksp,mpjac->submats,mpjac->submats,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
1320     }
1321   }
1322 
1323   if (pc->setfromoptionscalled){
1324     ierr = KSPSetFromOptions(mpjac->ksp);CHKERRQ(ierr);
1325   }
1326   PetscFunctionReturn(0);
1327 }
1328