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