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