xref: /petsc/src/ksp/pc/impls/bjacobi/bjacobi.c (revision c2efdce3a8b4106b7e14f0dd7333a5b437b351ad)
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   PC_BJacobi_Singleblock *bjac;
837   PetscBool              wasSetup;
838 
839   PetscFunctionBegin;
840 
841   if (!pc->setupcalled) {
842     const char *prefix;
843 
844     if (!jac->ksp) {
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       pc->ops->reset               = PCReset_BJacobi_Singleblock;
855       pc->ops->destroy             = PCDestroy_BJacobi_Singleblock;
856       pc->ops->apply               = PCApply_BJacobi_Singleblock;
857       pc->ops->applysymmetricleft  = PCApplySymmetricLeft_BJacobi_Singleblock;
858       pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Singleblock;
859       pc->ops->applytranspose      = PCApplyTranspose_BJacobi_Singleblock;
860       pc->ops->setuponblocks       = PCSetUpOnBlocks_BJacobi_Singleblock;
861 
862       ierr = PetscMalloc(sizeof(KSP),&jac->ksp);CHKERRQ(ierr);
863       jac->ksp[0] = ksp;
864 
865       ierr = PetscNewLog(pc,PC_BJacobi_Singleblock,&bjac);CHKERRQ(ierr);
866       jac->data    = (void*)bjac;
867     } else {
868       ksp  = jac->ksp[0];
869       bjac = (PC_BJacobi_Singleblock *)jac->data;
870     }
871 
872     /*
873       The reason we need to generate these vectors is to serve
874       as the right-hand side and solution vector for the solve on the
875       block. We do not need to allocate space for the vectors since
876       that is provided via VecPlaceArray() just before the call to
877       KSPSolve() on the block.
878     */
879     ierr = MatGetSize(pmat,&m,&m);CHKERRQ(ierr);
880     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,m,PETSC_NULL,&bjac->x);CHKERRQ(ierr);
881     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,m,PETSC_NULL,&bjac->y);CHKERRQ(ierr);
882     ierr = PetscLogObjectParent(pc,bjac->x);CHKERRQ(ierr);
883     ierr = PetscLogObjectParent(pc,bjac->y);CHKERRQ(ierr);
884   } else {
885     wasSetup = PETSC_TRUE;
886     ksp = jac->ksp[0];
887     bjac = (PC_BJacobi_Singleblock *)jac->data;
888   }
889   if (jac->use_true_local) {
890     ierr = KSPSetOperators(ksp,mat,pmat,pc->flag);CHKERRQ(ierr);
891   }  else {
892     ierr = KSPSetOperators(ksp,pmat,pmat,pc->flag);CHKERRQ(ierr);
893   }
894   if (!wasSetup && pc->setfromoptionscalled) {
895     ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
896   }
897   PetscFunctionReturn(0);
898 }
899 
900 /* ---------------------------------------------------------------------------------------------*/
901 #undef __FUNCT__
902 #define __FUNCT__ "PCReset_BJacobi_Multiblock"
903 PetscErrorCode PCReset_BJacobi_Multiblock(PC pc)
904 {
905   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
906   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
907   PetscErrorCode        ierr;
908   PetscInt              i;
909 
910   PetscFunctionBegin;
911   if (bjac && bjac->pmat) {
912     ierr = MatDestroyMatrices(jac->n_local,&bjac->pmat);CHKERRQ(ierr);
913     if (jac->use_true_local) {
914       ierr = MatDestroyMatrices(jac->n_local,&bjac->mat);CHKERRQ(ierr);
915     }
916   }
917 
918   for (i=0; i<jac->n_local; i++) {
919     ierr = KSPReset(jac->ksp[i]);CHKERRQ(ierr);
920     if (bjac && bjac->x) {
921       ierr = VecDestroy(&bjac->x[i]);CHKERRQ(ierr);
922       ierr = VecDestroy(&bjac->y[i]);CHKERRQ(ierr);
923       ierr = ISDestroy(&bjac->is[i]);CHKERRQ(ierr);
924     }
925   }
926   PetscFunctionReturn(0);
927 }
928 
929 #undef __FUNCT__
930 #define __FUNCT__ "PCDestroy_BJacobi_Multiblock"
931 PetscErrorCode PCDestroy_BJacobi_Multiblock(PC pc)
932 {
933   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
934   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
935   PetscErrorCode        ierr;
936   PetscInt              i;
937 
938   PetscFunctionBegin;
939   ierr = PCReset_BJacobi_Multiblock(pc);CHKERRQ(ierr);
940   if (bjac) {
941     ierr = PetscFree2(bjac->x,bjac->y);CHKERRQ(ierr);
942     ierr = PetscFree(bjac->starts);CHKERRQ(ierr);
943     ierr = PetscFree(bjac->is);CHKERRQ(ierr);
944   }
945   ierr = PetscFree(jac->data);CHKERRQ(ierr);
946   ierr = PetscFree(jac->l_lens);CHKERRQ(ierr);
947   ierr = PetscFree(jac->g_lens);CHKERRQ(ierr);
948   for (i=0; i<jac->n_local; i++) {
949     ierr = KSPDestroy(&jac->ksp[i]);CHKERRQ(ierr);
950   }
951   ierr = PetscFree(jac->ksp);CHKERRQ(ierr);
952   ierr = PetscFree(pc->data);CHKERRQ(ierr);
953   PetscFunctionReturn(0);
954 }
955 
956 #undef __FUNCT__
957 #define __FUNCT__ "PCSetUpOnBlocks_BJacobi_Multiblock"
958 PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiblock(PC pc)
959 {
960   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;
961   PetscErrorCode ierr;
962   PetscInt       i,n_local = jac->n_local;
963 
964   PetscFunctionBegin;
965   for (i=0; i<n_local; i++) {
966     ierr = KSPSetUp(jac->ksp[i]);CHKERRQ(ierr);
967   }
968   PetscFunctionReturn(0);
969 }
970 
971 /*
972       Preconditioner for block Jacobi
973 */
974 #undef __FUNCT__
975 #define __FUNCT__ "PCApply_BJacobi_Multiblock"
976 PetscErrorCode PCApply_BJacobi_Multiblock(PC pc,Vec x,Vec y)
977 {
978   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
979   PetscErrorCode        ierr;
980   PetscInt              i,n_local = jac->n_local;
981   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
982   PetscScalar           *xin,*yin;
983 
984   PetscFunctionBegin;
985   ierr = VecGetArray(x,&xin);CHKERRQ(ierr);
986   ierr = VecGetArray(y,&yin);CHKERRQ(ierr);
987   for (i=0; i<n_local; i++) {
988     /*
989        To avoid copying the subvector from x into a workspace we instead
990        make the workspace vector array point to the subpart of the array of
991        the global vector.
992     */
993     ierr = VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);CHKERRQ(ierr);
994     ierr = VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);CHKERRQ(ierr);
995 
996     ierr = PetscLogEventBegin(PC_SetUpOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);CHKERRQ(ierr);
997     ierr = KSPSolve(jac->ksp[i],bjac->x[i],bjac->y[i]);CHKERRQ(ierr);
998     ierr = PetscLogEventEnd(PC_SetUpOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);CHKERRQ(ierr);
999 
1000     ierr = VecResetArray(bjac->x[i]);CHKERRQ(ierr);
1001     ierr = VecResetArray(bjac->y[i]);CHKERRQ(ierr);
1002   }
1003   ierr = VecRestoreArray(x,&xin);CHKERRQ(ierr);
1004   ierr = VecRestoreArray(y,&yin);CHKERRQ(ierr);
1005   PetscFunctionReturn(0);
1006 }
1007 
1008 /*
1009       Preconditioner for block Jacobi
1010 */
1011 #undef __FUNCT__
1012 #define __FUNCT__ "PCApplyTranspose_BJacobi_Multiblock"
1013 PetscErrorCode PCApplyTranspose_BJacobi_Multiblock(PC pc,Vec x,Vec y)
1014 {
1015   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
1016   PetscErrorCode        ierr;
1017   PetscInt              i,n_local = jac->n_local;
1018   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
1019   PetscScalar           *xin,*yin;
1020 
1021   PetscFunctionBegin;
1022   ierr = VecGetArray(x,&xin);CHKERRQ(ierr);
1023   ierr = VecGetArray(y,&yin);CHKERRQ(ierr);
1024   for (i=0; i<n_local; i++) {
1025     /*
1026        To avoid copying the subvector from x into a workspace we instead
1027        make the workspace vector array point to the subpart of the array of
1028        the global vector.
1029     */
1030     ierr = VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);CHKERRQ(ierr);
1031     ierr = VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);CHKERRQ(ierr);
1032 
1033     ierr = PetscLogEventBegin(PC_ApplyTransposeOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);CHKERRQ(ierr);
1034     ierr = KSPSolveTranspose(jac->ksp[i],bjac->x[i],bjac->y[i]);CHKERRQ(ierr);
1035     ierr = PetscLogEventEnd(PC_ApplyTransposeOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);CHKERRQ(ierr);
1036   }
1037   ierr = VecRestoreArray(x,&xin);CHKERRQ(ierr);
1038   ierr = VecRestoreArray(y,&yin);CHKERRQ(ierr);
1039   PetscFunctionReturn(0);
1040 }
1041 
1042 #undef __FUNCT__
1043 #define __FUNCT__ "PCSetUp_BJacobi_Multiblock"
1044 static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC pc,Mat mat,Mat pmat)
1045 {
1046   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
1047   PetscErrorCode         ierr;
1048   PetscInt               m,n_local,N,M,start,i;
1049   const char             *prefix,*pprefix,*mprefix;
1050   KSP                    ksp;
1051   Vec                    x,y;
1052   PC_BJacobi_Multiblock  *bjac = (PC_BJacobi_Multiblock*)jac->data;
1053   PC                     subpc;
1054   IS                     is;
1055   MatReuse               scall;
1056 
1057   PetscFunctionBegin;
1058   ierr = MatGetLocalSize(pc->pmat,&M,&N);CHKERRQ(ierr);
1059 
1060   n_local = jac->n_local;
1061 
1062   if (jac->use_true_local) {
1063     PetscBool  same;
1064     ierr = PetscTypeCompare((PetscObject)mat,((PetscObject)pmat)->type_name,&same);CHKERRQ(ierr);
1065     if (!same) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Matrices not of same type");
1066   }
1067 
1068   if (!pc->setupcalled) {
1069     scall                  = MAT_INITIAL_MATRIX;
1070 
1071     if (!jac->ksp) {
1072       pc->ops->reset         = PCReset_BJacobi_Multiblock;
1073       pc->ops->destroy       = PCDestroy_BJacobi_Multiblock;
1074       pc->ops->apply         = PCApply_BJacobi_Multiblock;
1075       pc->ops->applytranspose= PCApplyTranspose_BJacobi_Multiblock;
1076       pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiblock;
1077 
1078       ierr = PetscNewLog(pc,PC_BJacobi_Multiblock,&bjac);CHKERRQ(ierr);
1079       ierr = PetscMalloc(n_local*sizeof(KSP),&jac->ksp);CHKERRQ(ierr);
1080       ierr = PetscLogObjectMemory(pc,sizeof(n_local*sizeof(KSP)));CHKERRQ(ierr);
1081       ierr = PetscMalloc2(n_local,Vec,&bjac->x,n_local,Vec,&bjac->y);CHKERRQ(ierr);
1082       ierr = PetscMalloc(n_local*sizeof(PetscScalar),&bjac->starts);CHKERRQ(ierr);
1083       ierr = PetscLogObjectMemory(pc,sizeof(n_local*sizeof(PetscScalar)));CHKERRQ(ierr);
1084 
1085       jac->data    = (void*)bjac;
1086       ierr = PetscMalloc(n_local*sizeof(IS),&bjac->is);CHKERRQ(ierr);
1087       ierr = PetscLogObjectMemory(pc,sizeof(n_local*sizeof(IS)));CHKERRQ(ierr);
1088 
1089       for (i=0; i<n_local; i++) {
1090         ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr);
1091         ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1092         ierr = PetscLogObjectParent(pc,ksp);CHKERRQ(ierr);
1093         ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
1094         ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr);
1095         ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
1096         ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
1097         ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
1098         jac->ksp[i]    = ksp;
1099       }
1100     } else {
1101       bjac = (PC_BJacobi_Multiblock*)jac->data;
1102     }
1103 
1104     start = 0;
1105     for (i=0; i<n_local; i++) {
1106         m = jac->l_lens[i];
1107         printf("m %d\n",m);
1108       /*
1109       The reason we need to generate these vectors is to serve
1110       as the right-hand side and solution vector for the solve on the
1111       block. We do not need to allocate space for the vectors since
1112       that is provided via VecPlaceArray() just before the call to
1113       KSPSolve() on the block.
1114 
1115       */
1116       ierr = VecCreateSeq(PETSC_COMM_SELF,m,&x);CHKERRQ(ierr);
1117       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,m,PETSC_NULL,&y);CHKERRQ(ierr);
1118       ierr = PetscLogObjectParent(pc,x);CHKERRQ(ierr);
1119       ierr = PetscLogObjectParent(pc,y);CHKERRQ(ierr);
1120       bjac->x[i]      = x;
1121       bjac->y[i]      = y;
1122       bjac->starts[i] = start;
1123 
1124       ierr = ISCreateStride(PETSC_COMM_SELF,m,start,1,&is);CHKERRQ(ierr);
1125       bjac->is[i] = is;
1126       ierr = PetscLogObjectParent(pc,is);CHKERRQ(ierr);
1127 
1128       start += m;
1129     }
1130   } else {
1131     bjac = (PC_BJacobi_Multiblock*)jac->data;
1132     /*
1133        Destroy the blocks from the previous iteration
1134     */
1135     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
1136       ierr = MatDestroyMatrices(n_local,&bjac->pmat);CHKERRQ(ierr);
1137       if (jac->use_true_local) {
1138         ierr = MatDestroyMatrices(n_local,&bjac->mat);CHKERRQ(ierr);
1139       }
1140       scall = MAT_INITIAL_MATRIX;
1141     } else {
1142       scall = MAT_REUSE_MATRIX;
1143     }
1144   }
1145 
1146   ierr = MatGetSubMatrices(pmat,n_local,bjac->is,bjac->is,scall,&bjac->pmat);CHKERRQ(ierr);
1147   if (jac->use_true_local) {
1148     ierr = PetscObjectGetOptionsPrefix((PetscObject)mat,&mprefix);CHKERRQ(ierr);
1149     ierr = MatGetSubMatrices(mat,n_local,bjac->is,bjac->is,scall,&bjac->mat);CHKERRQ(ierr);
1150   }
1151   /* Return control to the user so that the submatrices can be modified (e.g., to apply
1152      different boundary conditions for the submatrices than for the global problem) */
1153   ierr = PCModifySubMatrices(pc,n_local,bjac->is,bjac->is,bjac->pmat,pc->modifysubmatricesP);CHKERRQ(ierr);
1154 
1155   ierr = PetscObjectGetOptionsPrefix((PetscObject)pmat,&pprefix);CHKERRQ(ierr);
1156   for (i=0; i<n_local; i++) {
1157     ierr = PetscLogObjectParent(pc,bjac->pmat[i]);CHKERRQ(ierr);
1158     ierr = PetscObjectSetOptionsPrefix((PetscObject)bjac->pmat[i],pprefix);CHKERRQ(ierr);
1159     if (jac->use_true_local) {
1160       ierr = PetscLogObjectParent(pc,bjac->mat[i]);CHKERRQ(ierr);
1161       ierr = PetscObjectSetOptionsPrefix((PetscObject)bjac->mat[i],mprefix);CHKERRQ(ierr);
1162       ierr = KSPSetOperators(jac->ksp[i],bjac->mat[i],bjac->pmat[i],pc->flag);CHKERRQ(ierr);
1163     } else {
1164       ierr = KSPSetOperators(jac->ksp[i],bjac->pmat[i],bjac->pmat[i],pc->flag);CHKERRQ(ierr);
1165     }
1166     if(pc->setfromoptionscalled){
1167       ierr = KSPSetFromOptions(jac->ksp[i]);CHKERRQ(ierr);
1168     }
1169   }
1170   PetscFunctionReturn(0);
1171 }
1172 
1173 /* ---------------------------------------------------------------------------------------------*/
1174 /*
1175       These are for a single block with multiple processes;
1176 */
1177 #undef __FUNCT__
1178 #define __FUNCT__ "PCReset_BJacobi_Multiproc"
1179 static PetscErrorCode PCReset_BJacobi_Multiproc(PC pc)
1180 {
1181   PC_BJacobi           *jac = (PC_BJacobi*)pc->data;
1182   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1183   PetscErrorCode       ierr;
1184 
1185   PetscFunctionBegin;
1186   ierr = VecDestroy(&mpjac->ysub);CHKERRQ(ierr);
1187   ierr = VecDestroy(&mpjac->xsub);CHKERRQ(ierr);
1188   ierr = MatDestroy(&mpjac->submats);CHKERRQ(ierr);
1189   if (mpjac->ksp){ierr = KSPReset(mpjac->ksp);CHKERRQ(ierr);}
1190   PetscFunctionReturn(0);
1191 }
1192 
1193 #undef __FUNCT__
1194 #define __FUNCT__ "PCDestroy_BJacobi_Multiproc"
1195 static PetscErrorCode PCDestroy_BJacobi_Multiproc(PC pc)
1196 {
1197   PC_BJacobi           *jac = (PC_BJacobi*)pc->data;
1198   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1199   PetscErrorCode       ierr;
1200 
1201   PetscFunctionBegin;
1202   ierr = PCReset_BJacobi_Multiproc(pc);CHKERRQ(ierr);
1203   ierr = KSPDestroy(&mpjac->ksp);CHKERRQ(ierr);
1204   ierr = PetscSubcommDestroy(&mpjac->psubcomm);CHKERRQ(ierr);
1205 
1206   ierr = PetscFree(mpjac);CHKERRQ(ierr);
1207   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1208   PetscFunctionReturn(0);
1209 }
1210 
1211 #undef __FUNCT__
1212 #define __FUNCT__ "PCApply_BJacobi_Multiproc"
1213 static PetscErrorCode PCApply_BJacobi_Multiproc(PC pc,Vec x,Vec y)
1214 {
1215   PC_BJacobi           *jac = (PC_BJacobi*)pc->data;
1216   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1217   PetscErrorCode       ierr;
1218   PetscScalar          *xarray,*yarray;
1219 
1220   PetscFunctionBegin;
1221   /* place x's and y's local arrays into xsub and ysub */
1222   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
1223   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
1224   ierr = VecPlaceArray(mpjac->xsub,xarray);CHKERRQ(ierr);
1225   ierr = VecPlaceArray(mpjac->ysub,yarray);CHKERRQ(ierr);
1226 
1227   /* apply preconditioner on each matrix block */
1228   ierr = KSPSolve(mpjac->ksp,mpjac->xsub,mpjac->ysub);CHKERRQ(ierr);
1229 
1230   ierr = VecResetArray(mpjac->xsub);CHKERRQ(ierr);
1231   ierr = VecResetArray(mpjac->ysub);CHKERRQ(ierr);
1232   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
1233   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
1234   PetscFunctionReturn(0);
1235 }
1236 
1237 extern PetscErrorCode MatGetMultiProcBlock_MPIAIJ(Mat,MPI_Comm,Mat*);
1238 #undef __FUNCT__
1239 #define __FUNCT__ "PCSetUp_BJacobi_Multiproc"
1240 static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC pc)
1241 {
1242   PC_BJacobi           *jac = (PC_BJacobi*)pc->data;
1243   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1244   PetscErrorCode       ierr;
1245   PetscInt             m,n;
1246   MPI_Comm             comm = ((PetscObject)pc)->comm,subcomm=0;
1247   const char           *prefix;
1248 
1249   PetscFunctionBegin;
1250   if (jac->n_local > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only a single block in a subcommunicator is supported");
1251   jac->n_local = 1; /* currently only a single block is supported for a subcommunicator */
1252   if (!pc->setupcalled) {
1253     ierr = PetscNewLog(pc,PC_BJacobi_Multiproc,&mpjac);CHKERRQ(ierr);
1254     jac->data = (void*)mpjac;
1255 
1256     /* initialize datastructure mpjac */
1257     if (!jac->psubcomm) {
1258       /* Create default contiguous subcommunicatiors if user does not provide them */
1259       ierr = PetscSubcommCreate(comm,&jac->psubcomm);CHKERRQ(ierr);
1260       ierr = PetscSubcommSetNumber(jac->psubcomm,jac->n);CHKERRQ(ierr);
1261       ierr = PetscSubcommSetType(jac->psubcomm,PETSC_SUBCOMM_CONTIGUOUS);CHKERRQ(ierr);
1262       ierr = PetscLogObjectMemory(pc,sizeof(PetscSubcomm));CHKERRQ(ierr);
1263     }
1264     mpjac->psubcomm = jac->psubcomm;
1265     subcomm         = mpjac->psubcomm->comm;
1266 
1267     /* Get matrix blocks of pmat */
1268     ierr = MatGetMultiProcBlock_MPIAIJ(pc->pmat,subcomm,&mpjac->submats);CHKERRQ(ierr);
1269 
1270     /* create a new PC that processors in each subcomm have copy of */
1271     ierr = KSPCreate(subcomm,&mpjac->ksp);CHKERRQ(ierr);
1272     ierr = PetscObjectIncrementTabLevel((PetscObject)mpjac->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1273     ierr = PetscLogObjectParent(pc,mpjac->ksp);CHKERRQ(ierr);
1274     ierr = KSPSetOperators(mpjac->ksp,mpjac->submats,mpjac->submats,pc->flag);CHKERRQ(ierr);
1275     ierr = KSPGetPC(mpjac->ksp,&mpjac->pc);CHKERRQ(ierr);
1276 
1277     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
1278     ierr = KSPSetOptionsPrefix(mpjac->ksp,prefix);CHKERRQ(ierr);
1279     ierr = KSPAppendOptionsPrefix(mpjac->ksp,"sub_");CHKERRQ(ierr);
1280     /*
1281       PetscMPIInt rank,subsize,subrank;
1282       ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1283       ierr = MPI_Comm_size(subcomm,&subsize);CHKERRQ(ierr);
1284       ierr = MPI_Comm_rank(subcomm,&subrank);CHKERRQ(ierr);
1285 
1286       ierr = MatGetLocalSize(mpjac->submats,&m,PETSC_NULL);CHKERRQ(ierr);
1287       ierr = MatGetSize(mpjac->submats,&n,PETSC_NULL);CHKERRQ(ierr);
1288       ierr = PetscSynchronizedPrintf(comm,"[%d], sub-size %d,sub-rank %d\n",rank,subsize,subrank);
1289       ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr);
1290     */
1291 
1292     /* create dummy vectors xsub and ysub */
1293     ierr = MatGetLocalSize(mpjac->submats,&m,&n);CHKERRQ(ierr);
1294     ierr = VecCreateMPIWithArray(subcomm,n,PETSC_DECIDE,PETSC_NULL,&mpjac->xsub);CHKERRQ(ierr);
1295     ierr = VecCreateMPIWithArray(subcomm,m,PETSC_DECIDE,PETSC_NULL,&mpjac->ysub);CHKERRQ(ierr);
1296     ierr = PetscLogObjectParent(pc,mpjac->xsub);CHKERRQ(ierr);
1297     ierr = PetscLogObjectParent(pc,mpjac->ysub);CHKERRQ(ierr);
1298 
1299     pc->ops->reset   = PCReset_BJacobi_Multiproc;
1300     pc->ops->destroy = PCDestroy_BJacobi_Multiproc;
1301     pc->ops->apply   = PCApply_BJacobi_Multiproc;
1302   }
1303 
1304   if (pc->setupcalled && pc->flag == DIFFERENT_NONZERO_PATTERN) {
1305     /* destroy old matrix blocks, then get new matrix blocks */
1306     if (mpjac->submats) {
1307       ierr = MatDestroy(&mpjac->submats);CHKERRQ(ierr);
1308       subcomm = mpjac->psubcomm->comm;
1309       ierr = MatGetMultiProcBlock_MPIAIJ(pc->pmat,subcomm,&mpjac->submats);CHKERRQ(ierr);
1310       ierr = KSPSetOperators(mpjac->ksp,mpjac->submats,mpjac->submats,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
1311     }
1312   }
1313 
1314   if (pc->setfromoptionscalled){
1315     ierr = KSPSetFromOptions(mpjac->ksp);CHKERRQ(ierr);
1316   }
1317   PetscFunctionReturn(0);
1318 }
1319