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