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