Actual source code: bjacobi.c

petsc-3.11.4 2019-09-28
Report Typos and Errors

  2: /*
  3:    Defines a block Jacobi preconditioner.
  4: */

  6:  #include <../src/ksp/pc/impls/bjacobi/bjacobi.h>

  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);

 12: static PetscErrorCode PCSetUp_BJacobi(PC pc)
 13: {
 14:   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;
 15:   Mat            mat  = pc->mat,pmat = pc->pmat;
 17:   PetscBool      hasop;
 18:   PetscInt       N,M,start,i,sum,end;
 19:   PetscInt       bs,i_start=-1,i_end=-1;
 20:   PetscMPIInt    rank,size;
 21:   const char     *pprefix,*mprefix;

 24:   MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);
 25:   MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
 26:   MatGetLocalSize(pc->pmat,&M,&N);
 27:   MatGetBlockSize(pc->pmat,&bs);

 29:   if (jac->n > 0 && jac->n < size) {
 30:     PCSetUp_BJacobi_Multiproc(pc);
 31:     return(0);
 32:   }

 34:   /* --------------------------------------------------------------------------
 35:       Determines the number of blocks assigned to each processor
 36:   -----------------------------------------------------------------------------*/

 38:   /*   local block count  given */
 39:   if (jac->n_local > 0 && jac->n < 0) {
 40:     MPIU_Allreduce(&jac->n_local,&jac->n,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));
 41:     if (jac->l_lens) { /* check that user set these correctly */
 42:       sum = 0;
 43:       for (i=0; i<jac->n_local; i++) {
 44:         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");
 45:         sum += jac->l_lens[i];
 46:       }
 47:       if (sum != M) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local lens set incorrectly");
 48:     } else {
 49:       PetscMalloc1(jac->n_local,&jac->l_lens);
 50:       for (i=0; i<jac->n_local; i++) jac->l_lens[i] = bs*((M/bs)/jac->n_local + (((M/bs) % jac->n_local) > i));
 51:     }
 52:   } else if (jac->n > 0 && jac->n_local < 0) { /* global block count given */
 53:     /* global blocks given: determine which ones are local */
 54:     if (jac->g_lens) {
 55:       /* check if the g_lens is has valid entries */
 56:       for (i=0; i<jac->n; i++) {
 57:         if (!jac->g_lens[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Zero block not allowed");
 58:         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");
 59:       }
 60:       if (size == 1) {
 61:         jac->n_local = jac->n;
 62:         PetscMalloc1(jac->n_local,&jac->l_lens);
 63:         PetscMemcpy(jac->l_lens,jac->g_lens,jac->n_local*sizeof(PetscInt));
 64:         /* check that user set these correctly */
 65:         sum = 0;
 66:         for (i=0; i<jac->n_local; i++) sum += jac->l_lens[i];
 67:         if (sum != M) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Global lens set incorrectly");
 68:       } else {
 69:         MatGetOwnershipRange(pc->pmat,&start,&end);
 70:         /* loop over blocks determing first one owned by me */
 71:         sum = 0;
 72:         for (i=0; i<jac->n+1; i++) {
 73:           if (sum == start) { i_start = i; goto start_1;}
 74:           if (i < jac->n) sum += jac->g_lens[i];
 75:         }
 76:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Block sizes used in PCBJacobiSetTotalBlocks()\nare not compatible with parallel matrix layout");
 77: start_1:
 78:         for (i=i_start; i<jac->n+1; i++) {
 79:           if (sum == end) { i_end = i; goto end_1; }
 80:           if (i < jac->n) sum += jac->g_lens[i];
 81:         }
 82:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Block sizes used in PCBJacobiSetTotalBlocks()\nare not compatible with parallel matrix layout");
 83: end_1:
 84:         jac->n_local = i_end - i_start;
 85:         PetscMalloc1(jac->n_local,&jac->l_lens);
 86:         PetscMemcpy(jac->l_lens,jac->g_lens+i_start,jac->n_local*sizeof(PetscInt));
 87:       }
 88:     } else { /* no global blocks given, determine then using default layout */
 89:       jac->n_local = jac->n/size + ((jac->n % size) > rank);
 90:       PetscMalloc1(jac->n_local,&jac->l_lens);
 91:       for (i=0; i<jac->n_local; i++) {
 92:         jac->l_lens[i] = ((M/bs)/jac->n_local + (((M/bs) % jac->n_local) > i))*bs;
 93:         if (!jac->l_lens[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Too many blocks given");
 94:       }
 95:     }
 96:   } else if (jac->n < 0 && jac->n_local < 0) { /* no blocks given */
 97:     jac->n         = size;
 98:     jac->n_local   = 1;
 99:     PetscMalloc1(1,&jac->l_lens);
100:     jac->l_lens[0] = M;
101:   } else { /* jac->n > 0 && jac->n_local > 0 */
102:     if (!jac->l_lens) {
103:       PetscMalloc1(jac->n_local,&jac->l_lens);
104:       for (i=0; i<jac->n_local; i++) jac->l_lens[i] = bs*((M/bs)/jac->n_local + (((M/bs) % jac->n_local) > i));
105:     }
106:   }
107:   if (jac->n_local < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of blocks is less than number of processors");

109:   /* -------------------------
110:       Determines mat and pmat
111:   ---------------------------*/
112:   MatHasOperation(pc->mat,MATOP_GET_DIAGONAL_BLOCK,&hasop);
113:   if (!hasop && size == 1) {
114:     mat  = pc->mat;
115:     pmat = pc->pmat;
116:   } else {
117:     if (pc->useAmat) {
118:       /* use block from Amat matrix, not Pmat for local MatMult() */
119:       MatGetDiagonalBlock(pc->mat,&mat);
120:       /* make submatrix have same prefix as entire matrix */
121:       PetscObjectGetOptionsPrefix((PetscObject)pc->mat,&mprefix);
122:       PetscObjectSetOptionsPrefix((PetscObject)mat,mprefix);
123:     }
124:     if (pc->pmat != pc->mat || !pc->useAmat) {
125:       MatGetDiagonalBlock(pc->pmat,&pmat);
126:       /* make submatrix have same prefix as entire matrix */
127:       PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);
128:       PetscObjectSetOptionsPrefix((PetscObject)pmat,pprefix);
129:     } else pmat = mat;
130:   }

132:   /* ------
133:      Setup code depends on the number of blocks
134:   */
135:   if (jac->n_local == 1) {
136:     PCSetUp_BJacobi_Singleblock(pc,mat,pmat);
137:   } else {
138:     PCSetUp_BJacobi_Multiblock(pc,mat,pmat);
139:   }
140:   return(0);
141: }

143: /* Default destroy, if it has never been setup */
144: static PetscErrorCode PCDestroy_BJacobi(PC pc)
145: {
146:   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;

150:   PetscFree(jac->g_lens);
151:   PetscFree(jac->l_lens);
152:   PetscFree(pc->data);
153:   return(0);
154: }


157: static PetscErrorCode PCSetFromOptions_BJacobi(PetscOptionItems *PetscOptionsObject,PC pc)
158: {
159:   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;
161:   PetscInt       blocks,i;
162:   PetscBool      flg;

165:   PetscOptionsHead(PetscOptionsObject,"Block Jacobi options");
166:   PetscOptionsInt("-pc_bjacobi_blocks","Total number of blocks","PCBJacobiSetTotalBlocks",jac->n,&blocks,&flg);
167:   if (flg) {PCBJacobiSetTotalBlocks(pc,blocks,NULL);}
168:   PetscOptionsInt("-pc_bjacobi_local_blocks","Local number of blocks","PCBJacobiSetLocalBlocks",jac->n_local,&blocks,&flg);
169:   if (flg) {PCBJacobiSetLocalBlocks(pc,blocks,NULL);}
170:   if (jac->ksp) {
171:     /* The sub-KSP has already been set up (e.g., PCSetUp_BJacobi_Singleblock), but KSPSetFromOptions was not called
172:      * unless we had already been called. */
173:     for (i=0; i<jac->n_local; i++) {
174:       KSPSetFromOptions(jac->ksp[i]);
175:     }
176:   }
177:   PetscOptionsTail();
178:   return(0);
179: }

181:  #include <petscdraw.h>
182: static PetscErrorCode PCView_BJacobi(PC pc,PetscViewer viewer)
183: {
184:   PC_BJacobi           *jac   = (PC_BJacobi*)pc->data;
185:   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
186:   PetscErrorCode       ierr;
187:   PetscMPIInt          rank;
188:   PetscInt             i;
189:   PetscBool            iascii,isstring,isdraw;
190:   PetscViewer          sviewer;

193:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
194:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
195:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
196:   if (iascii) {
197:     if (pc->useAmat) {
198:       PetscViewerASCIIPrintf(viewer,"  using Amat local matrix, number of blocks = %D\n",jac->n);
199:     }
200:     PetscViewerASCIIPrintf(viewer,"  number of blocks = %D\n",jac->n);
201:     MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);
202:     if (jac->same_local_solves) {
203:       PetscViewerASCIIPrintf(viewer,"  Local solve is same for all blocks, in the following KSP and PC objects:\n");
204:       if (jac->ksp && !jac->psubcomm) {
205:         PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
206:         if (!rank) {
207:           PetscViewerASCIIPushTab(viewer);
208:           KSPView(jac->ksp[0],sviewer);
209:           PetscViewerASCIIPopTab(viewer);
210:         }
211:         PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
212:       } else if (mpjac && jac->ksp && mpjac->psubcomm) {
213:         PetscViewerGetSubViewer(viewer,mpjac->psubcomm->child,&sviewer);
214:         if (!mpjac->psubcomm->color) {
215:           PetscViewerASCIIPushTab(viewer);
216:           KSPView(*(jac->ksp),sviewer);
217:           PetscViewerASCIIPopTab(viewer);
218:         }
219:         PetscViewerRestoreSubViewer(viewer,mpjac->psubcomm->child,&sviewer);
220:       }
221:     } else {
222:       PetscInt n_global;
223:       MPIU_Allreduce(&jac->n_local,&n_global,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)pc));
224:       PetscViewerASCIIPushSynchronized(viewer);
225:       PetscViewerASCIIPrintf(viewer,"  Local solve info for each block is in the following KSP and PC objects:\n");
226:       PetscViewerASCIISynchronizedPrintf(viewer,"[%d] number of local blocks = %D, first local block number = %D\n",
227:                                                 rank,jac->n_local,jac->first_local);
228:       PetscViewerASCIIPushTab(viewer);
229:       PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
230:       for (i=0; i<jac->n_local; i++) {
231:         PetscViewerASCIISynchronizedPrintf(viewer,"[%d] local block number %D\n",rank,i);
232:         KSPView(jac->ksp[i],sviewer);
233:         PetscViewerASCIISynchronizedPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");
234:       }
235:       PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
236:       PetscViewerASCIIPopTab(viewer);
237:       PetscViewerFlush(viewer);
238:       PetscViewerASCIIPopSynchronized(viewer);
239:     }
240:   } else if (isstring) {
241:     PetscViewerStringSPrintf(viewer," blks=%D",jac->n);
242:     PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
243:     if (jac->ksp) {KSPView(jac->ksp[0],sviewer);}
244:     PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
245:   } else if (isdraw) {
246:     PetscDraw draw;
247:     char      str[25];
248:     PetscReal x,y,bottom,h;

250:     PetscViewerDrawGetDraw(viewer,0,&draw);
251:     PetscDrawGetCurrentPoint(draw,&x,&y);
252:     PetscSNPrintf(str,25,"Number blocks %D",jac->n);
253:     PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);
254:     bottom = y - h;
255:     PetscDrawPushCurrentPoint(draw,x,bottom);
256:     /* warning the communicator on viewer is different then on ksp in parallel */
257:     if (jac->ksp) {KSPView(jac->ksp[0],viewer);}
258:     PetscDrawPopCurrentPoint(draw);
259:   }
260:   return(0);
261: }

263: /* -------------------------------------------------------------------------------------*/

265: static PetscErrorCode  PCBJacobiGetSubKSP_BJacobi(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp)
266: {
267:   PC_BJacobi *jac = (PC_BJacobi*)pc->data;;

270:   if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call KSPSetUp() or PCSetUp() first");

272:   if (n_local) *n_local = jac->n_local;
273:   if (first_local) *first_local = jac->first_local;
274:   *ksp                   = jac->ksp;
275:   jac->same_local_solves = PETSC_FALSE;        /* Assume that local solves are now different;
276:                                                   not necessarily true though!  This flag is
277:                                                   used only for PCView_BJacobi() */
278:   return(0);
279: }

281: static PetscErrorCode  PCBJacobiSetTotalBlocks_BJacobi(PC pc,PetscInt blocks,PetscInt *lens)
282: {
283:   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;

287:   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");
288:   jac->n = blocks;
289:   if (!lens) jac->g_lens = 0;
290:   else {
291:     PetscMalloc1(blocks,&jac->g_lens);
292:     PetscLogObjectMemory((PetscObject)pc,blocks*sizeof(PetscInt));
293:     PetscMemcpy(jac->g_lens,lens,blocks*sizeof(PetscInt));
294:   }
295:   return(0);
296: }

298: static PetscErrorCode  PCBJacobiGetTotalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[])
299: {
300:   PC_BJacobi *jac = (PC_BJacobi*) pc->data;

303:   *blocks = jac->n;
304:   if (lens) *lens = jac->g_lens;
305:   return(0);
306: }

308: static PetscErrorCode  PCBJacobiSetLocalBlocks_BJacobi(PC pc,PetscInt blocks,const PetscInt lens[])
309: {
310:   PC_BJacobi     *jac;

314:   jac = (PC_BJacobi*)pc->data;

316:   jac->n_local = blocks;
317:   if (!lens) jac->l_lens = 0;
318:   else {
319:     PetscMalloc1(blocks,&jac->l_lens);
320:     PetscLogObjectMemory((PetscObject)pc,blocks*sizeof(PetscInt));
321:     PetscMemcpy(jac->l_lens,lens,blocks*sizeof(PetscInt));
322:   }
323:   return(0);
324: }

326: static PetscErrorCode  PCBJacobiGetLocalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[])
327: {
328:   PC_BJacobi *jac = (PC_BJacobi*) pc->data;

331:   *blocks = jac->n_local;
332:   if (lens) *lens = jac->l_lens;
333:   return(0);
334: }

336: /* -------------------------------------------------------------------------------------*/

338: /*@C
339:    PCBJacobiGetSubKSP - Gets the local KSP contexts for all blocks on
340:    this processor.

342:    Not Collective

344:    Input Parameter:
345: .  pc - the preconditioner context

347:    Output Parameters:
348: +  n_local - the number of blocks on this processor, or NULL
349: .  first_local - the global number of the first block on this processor, or NULL
350: -  ksp - the array of KSP contexts

352:    Notes:
353:    After PCBJacobiGetSubKSP() the array of KSP contexts is not to be freed.

355:    Currently for some matrix implementations only 1 block per processor
356:    is supported.

358:    You must call KSPSetUp() or PCSetUp() before calling PCBJacobiGetSubKSP().

360:    Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs.
361:       You can call PCBJacobiGetSubKSP(pc,nlocal,firstlocal,PETSC_NULL_KSP,ierr) to determine how large the
362:       KSP array must be.

364:    Level: advanced

366: .keywords:  block, Jacobi, get, sub, KSP, context

368: .seealso: PCBJacobiGetSubKSP()
369: @*/
370: PetscErrorCode  PCBJacobiGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
371: {

376:   PetscUseMethod(pc,"PCBJacobiGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));
377:   return(0);
378: }

380: /*@
381:    PCBJacobiSetTotalBlocks - Sets the global number of blocks for the block
382:    Jacobi preconditioner.

384:    Collective on PC

386:    Input Parameters:
387: +  pc - the preconditioner context
388: .  blocks - the number of blocks
389: -  lens - [optional] integer array containing the size of each block

391:    Options Database Key:
392: .  -pc_bjacobi_blocks <blocks> - Sets the number of global blocks

394:    Notes:
395:    Currently only a limited number of blocking configurations are supported.
396:    All processors sharing the PC must call this routine with the same data.

398:    Level: intermediate

400: .keywords:  set, number, Jacobi, global, total, blocks

402: .seealso: PCSetUseAmat(), PCBJacobiSetLocalBlocks()
403: @*/
404: PetscErrorCode  PCBJacobiSetTotalBlocks(PC pc,PetscInt blocks,const PetscInt lens[])
405: {

410:   if (blocks <= 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Must have positive blocks");
411:   PetscTryMethod(pc,"PCBJacobiSetTotalBlocks_C",(PC,PetscInt,const PetscInt[]),(pc,blocks,lens));
412:   return(0);
413: }

415: /*@C
416:    PCBJacobiGetTotalBlocks - Gets the global number of blocks for the block
417:    Jacobi preconditioner.

419:    Not Collective

421:    Input Parameter:
422: .  pc - the preconditioner context

424:    Output parameters:
425: +  blocks - the number of blocks
426: -  lens - integer array containing the size of each block

428:    Level: intermediate

430: .keywords:  get, number, Jacobi, global, total, blocks

432: .seealso: PCSetUseAmat(), PCBJacobiGetLocalBlocks()
433: @*/
434: PetscErrorCode  PCBJacobiGetTotalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
435: {

441:   PetscUseMethod(pc,"PCBJacobiGetTotalBlocks_C",(PC,PetscInt*, const PetscInt *[]),(pc,blocks,lens));
442:   return(0);
443: }

445: /*@
446:    PCBJacobiSetLocalBlocks - Sets the local number of blocks for the block
447:    Jacobi preconditioner.

449:    Not Collective

451:    Input Parameters:
452: +  pc - the preconditioner context
453: .  blocks - the number of blocks
454: -  lens - [optional] integer array containing size of each block

456:    Options Database Key:
457: .  -pc_bjacobi_local_blocks <blocks> - Sets the number of local blocks

459:    Note:
460:    Currently only a limited number of blocking configurations are supported.

462:    Level: intermediate

464: .keywords: PC, set, number, Jacobi, local, blocks

466: .seealso: PCSetUseAmat(), PCBJacobiSetTotalBlocks()
467: @*/
468: PetscErrorCode  PCBJacobiSetLocalBlocks(PC pc,PetscInt blocks,const PetscInt lens[])
469: {

474:   if (blocks < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have nonegative blocks");
475:   PetscTryMethod(pc,"PCBJacobiSetLocalBlocks_C",(PC,PetscInt,const PetscInt []),(pc,blocks,lens));
476:   return(0);
477: }

479: /*@C
480:    PCBJacobiGetLocalBlocks - Gets the local number of blocks for the block
481:    Jacobi preconditioner.

483:    Not Collective

485:    Input Parameters:
486: +  pc - the preconditioner context
487: .  blocks - the number of blocks
488: -  lens - [optional] integer array containing size of each block

490:    Note:
491:    Currently only a limited number of blocking configurations are supported.

493:    Level: intermediate

495: .keywords: PC, get, number, Jacobi, local, blocks

497: .seealso: PCSetUseAmat(), PCBJacobiGetTotalBlocks()
498: @*/
499: PetscErrorCode  PCBJacobiGetLocalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
500: {

506:   PetscUseMethod(pc,"PCBJacobiGetLocalBlocks_C",(PC,PetscInt*, const PetscInt *[]),(pc,blocks,lens));
507:   return(0);
508: }

510: /* -----------------------------------------------------------------------------------*/

512: /*MC
513:    PCBJACOBI - Use block Jacobi preconditioning, each block is (approximately) solved with
514:            its own KSP object.

516:    Options Database Keys:
517: +  -pc_use_amat - use Amat to apply block of operator in inner Krylov method
518: -  -pc_bjacobi_blocks <n> - use n total blocks

520:    Notes:
521:     Each processor can have one or more blocks, or a single block can be shared by several processes. Defaults to one block per processor.

523:      To set options on the solvers for each block append -sub_ to all the KSP, KSP, and PC
524:         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly

526:      To set the options on the solvers separate for each block call PCBJacobiGetSubKSP()
527:          and set the options directly on the resulting KSP object (you can access its PC
528:          KSPGetPC())

530:      For GPU-based vectors (CUDA, ViennaCL) it is recommended to use exactly one block per MPI process for best
531:          performance.  Different block partitioning may lead to additional data transfers
532:          between host and GPU that lead to degraded performance.

534:      The options prefix for each block is sub_, for example -sub_pc_type lu.

536:      When multiple processes share a single block, each block encompasses exactly all the unknowns owned its set of processes.

538:    Level: beginner

540:    Concepts: block Jacobi

542: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
543:            PCASM, PCSetUseAmat(), PCGetUseAmat(), PCBJacobiGetSubKSP(), PCBJacobiSetTotalBlocks(),
544:            PCBJacobiSetLocalBlocks(), PCSetModifySubmatrices()
545: M*/

547: PETSC_EXTERN PetscErrorCode PCCreate_BJacobi(PC pc)
548: {
550:   PetscMPIInt    rank;
551:   PC_BJacobi     *jac;

554:   PetscNewLog(pc,&jac);
555:   MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);

557:   pc->ops->apply           = 0;
558:   pc->ops->applytranspose  = 0;
559:   pc->ops->setup           = PCSetUp_BJacobi;
560:   pc->ops->destroy         = PCDestroy_BJacobi;
561:   pc->ops->setfromoptions  = PCSetFromOptions_BJacobi;
562:   pc->ops->view            = PCView_BJacobi;
563:   pc->ops->applyrichardson = 0;

565:   pc->data               = (void*)jac;
566:   jac->n                 = -1;
567:   jac->n_local           = -1;
568:   jac->first_local       = rank;
569:   jac->ksp               = 0;
570:   jac->same_local_solves = PETSC_TRUE;
571:   jac->g_lens            = 0;
572:   jac->l_lens            = 0;
573:   jac->psubcomm          = 0;

575:   PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiGetSubKSP_C",PCBJacobiGetSubKSP_BJacobi);
576:   PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiSetTotalBlocks_C",PCBJacobiSetTotalBlocks_BJacobi);
577:   PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiGetTotalBlocks_C",PCBJacobiGetTotalBlocks_BJacobi);
578:   PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiSetLocalBlocks_C",PCBJacobiSetLocalBlocks_BJacobi);
579:   PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiGetLocalBlocks_C",PCBJacobiGetLocalBlocks_BJacobi);
580:   return(0);
581: }

583: /* --------------------------------------------------------------------------------------------*/
584: /*
585:         These are for a single block per processor; works for AIJ, BAIJ; Seq and MPI
586: */
587: static PetscErrorCode PCReset_BJacobi_Singleblock(PC pc)
588: {
589:   PC_BJacobi             *jac  = (PC_BJacobi*)pc->data;
590:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
591:   PetscErrorCode         ierr;

594:   KSPReset(jac->ksp[0]);
595:   VecDestroy(&bjac->x);
596:   VecDestroy(&bjac->y);
597:   return(0);
598: }

600: static PetscErrorCode PCDestroy_BJacobi_Singleblock(PC pc)
601: {
602:   PC_BJacobi             *jac  = (PC_BJacobi*)pc->data;
603:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
604:   PetscErrorCode         ierr;

607:   PCReset_BJacobi_Singleblock(pc);
608:   KSPDestroy(&jac->ksp[0]);
609:   PetscFree(jac->ksp);
610:   PetscFree(jac->l_lens);
611:   PetscFree(jac->g_lens);
612:   PetscFree(bjac);
613:   PetscFree(pc->data);
614:   return(0);
615: }

617: static PetscErrorCode PCSetUpOnBlocks_BJacobi_Singleblock(PC pc)
618: {
619:   PetscErrorCode     ierr;
620:   PC_BJacobi         *jac = (PC_BJacobi*)pc->data;
621:   KSP                subksp = jac->ksp[0];
622:   KSPConvergedReason reason;

625:   KSPSetUp(subksp);
626:   KSPGetConvergedReason(subksp,&reason);
627:   if (reason == KSP_DIVERGED_PC_FAILED) {
628:     pc->failedreason = PC_SUBPC_ERROR;
629:   }
630:   return(0);
631: }

633: static PetscErrorCode PCApply_BJacobi_Singleblock(PC pc,Vec x,Vec y)
634: {
635:   PetscErrorCode         ierr;
636:   PC_BJacobi             *jac  = (PC_BJacobi*)pc->data;
637:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;

640:   VecGetLocalVectorRead(x, bjac->x);
641:   VecGetLocalVector(y, bjac->y);
642:  /* Since the inner KSP matrix may point directly to the diagonal block of an MPI matrix the inner
643:      matrix may change even if the outter KSP/PC has not updated the preconditioner, this will trigger a rebuild
644:      of the inner preconditioner automatically unless we pass down the outter preconditioners reuse flag.*/
645:   KSPSetReusePreconditioner(jac->ksp[0],pc->reusepreconditioner);
646:   KSPSolve(jac->ksp[0],bjac->x,bjac->y);
647:   KSPCheckSolve(jac->ksp[0],pc,bjac->y);
648:   VecRestoreLocalVectorRead(x, bjac->x);
649:   VecRestoreLocalVector(y, bjac->y);
650:   return(0);
651: }

653: static PetscErrorCode PCApplySymmetricLeft_BJacobi_Singleblock(PC pc,Vec x,Vec y)
654: {
655:   PetscErrorCode         ierr;
656:   PC_BJacobi             *jac  = (PC_BJacobi*)pc->data;
657:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
658:   PetscScalar            *y_array;
659:   const PetscScalar      *x_array;
660:   PC                     subpc;

663:   /*
664:       The VecPlaceArray() is to avoid having to copy the
665:     y vector into the bjac->x vector. The reason for
666:     the bjac->x vector is that we need a sequential vector
667:     for the sequential solve.
668:   */
669:   VecGetArrayRead(x,&x_array);
670:   VecGetArray(y,&y_array);
671:   VecPlaceArray(bjac->x,x_array);
672:   VecPlaceArray(bjac->y,y_array);
673:   /* apply the symmetric left portion of the inner PC operator */
674:   /* note this by-passes the inner KSP and its options completely */
675:   KSPGetPC(jac->ksp[0],&subpc);
676:   PCApplySymmetricLeft(subpc,bjac->x,bjac->y);
677:   VecResetArray(bjac->x);
678:   VecResetArray(bjac->y);
679:   VecRestoreArrayRead(x,&x_array);
680:   VecRestoreArray(y,&y_array);
681:   return(0);
682: }

684: static PetscErrorCode PCApplySymmetricRight_BJacobi_Singleblock(PC pc,Vec x,Vec y)
685: {
686:   PetscErrorCode         ierr;
687:   PC_BJacobi             *jac  = (PC_BJacobi*)pc->data;
688:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
689:   PetscScalar            *y_array;
690:   const PetscScalar      *x_array;
691:   PC                     subpc;

694:   /*
695:       The VecPlaceArray() is to avoid having to copy the
696:     y vector into the bjac->x vector. The reason for
697:     the bjac->x vector is that we need a sequential vector
698:     for the sequential solve.
699:   */
700:   VecGetArrayRead(x,&x_array);
701:   VecGetArray(y,&y_array);
702:   VecPlaceArray(bjac->x,x_array);
703:   VecPlaceArray(bjac->y,y_array);

705:   /* apply the symmetric right portion of the inner PC operator */
706:   /* note this by-passes the inner KSP and its options completely */

708:   KSPGetPC(jac->ksp[0],&subpc);
709:   PCApplySymmetricRight(subpc,bjac->x,bjac->y);

711:   VecRestoreArrayRead(x,&x_array);
712:   VecRestoreArray(y,&y_array);
713:   return(0);
714: }

716: static PetscErrorCode PCApplyTranspose_BJacobi_Singleblock(PC pc,Vec x,Vec y)
717: {
718:   PetscErrorCode         ierr;
719:   PC_BJacobi             *jac  = (PC_BJacobi*)pc->data;
720:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
721:   PetscScalar            *y_array;
722:   const PetscScalar      *x_array;

725:   /*
726:       The VecPlaceArray() is to avoid having to copy the
727:     y vector into the bjac->x vector. The reason for
728:     the bjac->x vector is that we need a sequential vector
729:     for the sequential solve.
730:   */
731:   VecGetArrayRead(x,&x_array);
732:   VecGetArray(y,&y_array);
733:   VecPlaceArray(bjac->x,x_array);
734:   VecPlaceArray(bjac->y,y_array);
735:   KSPSolveTranspose(jac->ksp[0],bjac->x,bjac->y);
736:   KSPCheckSolve(jac->ksp[0],pc,bjac->y);
737:   VecResetArray(bjac->x);
738:   VecResetArray(bjac->y);
739:   VecRestoreArrayRead(x,&x_array);
740:   VecRestoreArray(y,&y_array);
741:   return(0);
742: }

744: static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC pc,Mat mat,Mat pmat)
745: {
746:   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
747:   PetscErrorCode         ierr;
748:   PetscInt               m;
749:   KSP                    ksp;
750:   PC_BJacobi_Singleblock *bjac;
751:   PetscBool              wasSetup = PETSC_TRUE;
752: #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_VIENNACL)
753:   PetscBool              is_gpumatrix = PETSC_FALSE;
754: #endif

757:   if (!pc->setupcalled) {
758:     const char *prefix;

760:     if (!jac->ksp) {
761:       wasSetup = PETSC_FALSE;

763:       KSPCreate(PETSC_COMM_SELF,&ksp);
764:       KSPSetErrorIfNotConverged(ksp,pc->erroriffailure);
765:       PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);
766:       PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);
767:       KSPSetType(ksp,KSPPREONLY);
768:       PCGetOptionsPrefix(pc,&prefix);
769:       KSPSetOptionsPrefix(ksp,prefix);
770:       KSPAppendOptionsPrefix(ksp,"sub_");

772:       pc->ops->reset               = PCReset_BJacobi_Singleblock;
773:       pc->ops->destroy             = PCDestroy_BJacobi_Singleblock;
774:       pc->ops->apply               = PCApply_BJacobi_Singleblock;
775:       pc->ops->applysymmetricleft  = PCApplySymmetricLeft_BJacobi_Singleblock;
776:       pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Singleblock;
777:       pc->ops->applytranspose      = PCApplyTranspose_BJacobi_Singleblock;
778:       pc->ops->setuponblocks       = PCSetUpOnBlocks_BJacobi_Singleblock;

780:       PetscMalloc1(1,&jac->ksp);
781:       jac->ksp[0] = ksp;

783:       PetscNewLog(pc,&bjac);
784:       jac->data = (void*)bjac;
785:     } else {
786:       ksp  = jac->ksp[0];
787:       bjac = (PC_BJacobi_Singleblock*)jac->data;
788:     }

790:     /*
791:       The reason we need to generate these vectors is to serve
792:       as the right-hand side and solution vector for the solve on the
793:       block. We do not need to allocate space for the vectors since
794:       that is provided via VecPlaceArray() just before the call to
795:       KSPSolve() on the block.
796:     */
797:     MatGetSize(pmat,&m,&m);
798:     VecCreateSeqWithArray(PETSC_COMM_SELF,1,m,NULL,&bjac->x);
799:     VecCreateSeqWithArray(PETSC_COMM_SELF,1,m,NULL,&bjac->y);
800: #if defined(PETSC_HAVE_CUDA)
801:     PetscObjectTypeCompareAny((PetscObject)pmat,&is_gpumatrix,MATAIJCUSPARSE,MATSEQAIJCUSPARSE,MATMPIAIJCUSPARSE,"");
802:     if (is_gpumatrix) {
803:       VecSetType(bjac->x,VECCUDA);
804:       VecSetType(bjac->y,VECCUDA);
805:     }
806: #endif
807: #if defined(PETSC_HAVE_VIENNACL)
808:     PetscObjectTypeCompareAny((PetscObject)pmat,&is_gpumatrix,MATAIJVIENNACL,MATSEQAIJVIENNACL,MATMPIAIJVIENNACL,"");
809:     if (is_gpumatrix) {
810:       VecSetType(bjac->x,VECVIENNACL);
811:       VecSetType(bjac->y,VECVIENNACL);
812:     }
813: #endif
814:     PetscLogObjectParent((PetscObject)pc,(PetscObject)bjac->x);
815:     PetscLogObjectParent((PetscObject)pc,(PetscObject)bjac->y);
816:   } else {
817:     ksp  = jac->ksp[0];
818:     bjac = (PC_BJacobi_Singleblock*)jac->data;
819:   }
820:   if (pc->useAmat) {
821:     KSPSetOperators(ksp,mat,pmat);
822:   } else {
823:     KSPSetOperators(ksp,pmat,pmat);
824:   }
825:   if (!wasSetup && pc->setfromoptionscalled) {
826:     /* If PCSetFromOptions_BJacobi is called later, KSPSetFromOptions will be called at that time. */
827:     KSPSetFromOptions(ksp);
828:   }
829:   return(0);
830: }

832: /* ---------------------------------------------------------------------------------------------*/
833: static PetscErrorCode PCReset_BJacobi_Multiblock(PC pc)
834: {
835:   PC_BJacobi            *jac  = (PC_BJacobi*)pc->data;
836:   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
837:   PetscErrorCode        ierr;
838:   PetscInt              i;

841:   if (bjac && bjac->pmat) {
842:     MatDestroyMatrices(jac->n_local,&bjac->pmat);
843:     if (pc->useAmat) {
844:       MatDestroyMatrices(jac->n_local,&bjac->mat);
845:     }
846:   }

848:   for (i=0; i<jac->n_local; i++) {
849:     KSPReset(jac->ksp[i]);
850:     if (bjac && bjac->x) {
851:       VecDestroy(&bjac->x[i]);
852:       VecDestroy(&bjac->y[i]);
853:       ISDestroy(&bjac->is[i]);
854:     }
855:   }
856:   PetscFree(jac->l_lens);
857:   PetscFree(jac->g_lens);
858:   return(0);
859: }

861: static PetscErrorCode PCDestroy_BJacobi_Multiblock(PC pc)
862: {
863:   PC_BJacobi            *jac  = (PC_BJacobi*)pc->data;
864:   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
865:   PetscErrorCode        ierr;
866:   PetscInt              i;

869:   PCReset_BJacobi_Multiblock(pc);
870:   if (bjac) {
871:     PetscFree2(bjac->x,bjac->y);
872:     PetscFree(bjac->starts);
873:     PetscFree(bjac->is);
874:   }
875:   PetscFree(jac->data);
876:   for (i=0; i<jac->n_local; i++) {
877:     KSPDestroy(&jac->ksp[i]);
878:   }
879:   PetscFree(jac->ksp);
880:   PetscFree(pc->data);
881:   return(0);
882: }

884: static PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiblock(PC pc)
885: {
886:   PC_BJacobi         *jac = (PC_BJacobi*)pc->data;
887:   PetscErrorCode     ierr;
888:   PetscInt           i,n_local = jac->n_local;
889:   KSPConvergedReason reason;

892:   for (i=0; i<n_local; i++) {
893:     KSPSetUp(jac->ksp[i]);
894:     KSPGetConvergedReason(jac->ksp[i],&reason);
895:     if (reason == KSP_DIVERGED_PC_FAILED) {
896:       pc->failedreason = PC_SUBPC_ERROR;
897:     }
898:   }
899:   return(0);
900: }

902: /*
903:       Preconditioner for block Jacobi
904: */
905: static PetscErrorCode PCApply_BJacobi_Multiblock(PC pc,Vec x,Vec y)
906: {
907:   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
908:   PetscErrorCode        ierr;
909:   PetscInt              i,n_local = jac->n_local;
910:   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
911:   PetscScalar           *yin;
912:   const PetscScalar     *xin;

915:   VecGetArrayRead(x,&xin);
916:   VecGetArray(y,&yin);
917:   for (i=0; i<n_local; i++) {
918:     /*
919:        To avoid copying the subvector from x into a workspace we instead
920:        make the workspace vector array point to the subpart of the array of
921:        the global vector.
922:     */
923:     VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);
924:     VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);

926:     PetscLogEventBegin(PC_ApplyOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);
927:     KSPSolve(jac->ksp[i],bjac->x[i],bjac->y[i]);
928:     KSPCheckSolve(jac->ksp[i],pc,bjac->y[i]);
929:     PetscLogEventEnd(PC_ApplyOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);

931:     VecResetArray(bjac->x[i]);
932:     VecResetArray(bjac->y[i]);
933:   }
934:   VecRestoreArrayRead(x,&xin);
935:   VecRestoreArray(y,&yin);
936:   return(0);
937: }

939: /*
940:       Preconditioner for block Jacobi
941: */
942: static PetscErrorCode PCApplyTranspose_BJacobi_Multiblock(PC pc,Vec x,Vec y)
943: {
944:   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
945:   PetscErrorCode        ierr;
946:   PetscInt              i,n_local = jac->n_local;
947:   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
948:   PetscScalar           *yin;
949:   const PetscScalar     *xin;

952:   VecGetArrayRead(x,&xin);
953:   VecGetArray(y,&yin);
954:   for (i=0; i<n_local; i++) {
955:     /*
956:        To avoid copying the subvector from x into a workspace we instead
957:        make the workspace vector array point to the subpart of the array of
958:        the global vector.
959:     */
960:     VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);
961:     VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);

963:     PetscLogEventBegin(PC_ApplyTransposeOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);
964:     KSPSolveTranspose(jac->ksp[i],bjac->x[i],bjac->y[i]);
965:     KSPCheckSolve(jac->ksp[i],pc,bjac->y[i]);
966:     PetscLogEventEnd(PC_ApplyTransposeOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);

968:     VecResetArray(bjac->x[i]);
969:     VecResetArray(bjac->y[i]);
970:   }
971:   VecRestoreArrayRead(x,&xin);
972:   VecRestoreArray(y,&yin);
973:   return(0);
974: }

976: static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC pc,Mat mat,Mat pmat)
977: {
978:   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
979:   PetscErrorCode        ierr;
980:   PetscInt              m,n_local,N,M,start,i;
981:   const char            *prefix,*pprefix,*mprefix;
982:   KSP                   ksp;
983:   Vec                   x,y;
984:   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
985:   PC                    subpc;
986:   IS                    is;
987:   MatReuse              scall;
988: #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_VIENNACL)
989:   PetscBool              is_gpumatrix = PETSC_FALSE;
990: #endif

993:   MatGetLocalSize(pc->pmat,&M,&N);

995:   n_local = jac->n_local;

997:   if (pc->useAmat) {
998:     PetscBool same;
999:     PetscObjectTypeCompare((PetscObject)mat,((PetscObject)pmat)->type_name,&same);
1000:     if (!same) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Matrices not of same type");
1001:   }

1003:   if (!pc->setupcalled) {
1004:     scall = MAT_INITIAL_MATRIX;

1006:     if (!jac->ksp) {
1007:       pc->ops->reset         = PCReset_BJacobi_Multiblock;
1008:       pc->ops->destroy       = PCDestroy_BJacobi_Multiblock;
1009:       pc->ops->apply         = PCApply_BJacobi_Multiblock;
1010:       pc->ops->applytranspose= PCApplyTranspose_BJacobi_Multiblock;
1011:       pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiblock;

1013:       PetscNewLog(pc,&bjac);
1014:       PetscMalloc1(n_local,&jac->ksp);
1015:       PetscLogObjectMemory((PetscObject)pc,sizeof(n_local*sizeof(KSP)));
1016:       PetscMalloc2(n_local,&bjac->x,n_local,&bjac->y);
1017:       PetscMalloc1(n_local,&bjac->starts);
1018:       PetscLogObjectMemory((PetscObject)pc,sizeof(n_local*sizeof(PetscScalar)));

1020:       jac->data = (void*)bjac;
1021:       PetscMalloc1(n_local,&bjac->is);
1022:       PetscLogObjectMemory((PetscObject)pc,sizeof(n_local*sizeof(IS)));

1024:       for (i=0; i<n_local; i++) {
1025:         KSPCreate(PETSC_COMM_SELF,&ksp);
1026:         KSPSetErrorIfNotConverged(ksp,pc->erroriffailure);
1027:         PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);
1028:         PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);
1029:         KSPSetType(ksp,KSPPREONLY);
1030:         KSPGetPC(ksp,&subpc);
1031:         PCGetOptionsPrefix(pc,&prefix);
1032:         KSPSetOptionsPrefix(ksp,prefix);
1033:         KSPAppendOptionsPrefix(ksp,"sub_");

1035:         jac->ksp[i] = ksp;
1036:       }
1037:     } else {
1038:       bjac = (PC_BJacobi_Multiblock*)jac->data;
1039:     }

1041:     start = 0;
1042:     for (i=0; i<n_local; i++) {
1043:       m = jac->l_lens[i];
1044:       /*
1045:       The reason we need to generate these vectors is to serve
1046:       as the right-hand side and solution vector for the solve on the
1047:       block. We do not need to allocate space for the vectors since
1048:       that is provided via VecPlaceArray() just before the call to
1049:       KSPSolve() on the block.

1051:       */
1052:       VecCreateSeq(PETSC_COMM_SELF,m,&x);
1053:       VecCreateSeqWithArray(PETSC_COMM_SELF,1,m,NULL,&y);
1054: #if defined(PETSC_HAVE_CUDA)
1055:       PetscObjectTypeCompareAny((PetscObject)pmat,&is_gpumatrix,MATAIJCUSPARSE,MATSEQAIJCUSPARSE,MATMPIAIJCUSPARSE,"");
1056:       if (is_gpumatrix) {
1057:         VecSetType(x,VECCUDA);
1058:         VecSetType(y,VECCUDA);
1059:       }
1060: #endif
1061: #if defined(PETSC_HAVE_VIENNACL)
1062:       PetscObjectTypeCompareAny((PetscObject)pmat,&is_gpumatrix,MATAIJVIENNACL,MATSEQAIJVIENNACL,MATMPIAIJVIENNACL,"");
1063:       if (is_gpumatrix) {
1064:         VecSetType(x,VECVIENNACL);
1065:         VecSetType(y,VECVIENNACL);
1066:       }
1067: #endif
1068:       PetscLogObjectParent((PetscObject)pc,(PetscObject)x);
1069:       PetscLogObjectParent((PetscObject)pc,(PetscObject)y);

1071:       bjac->x[i]      = x;
1072:       bjac->y[i]      = y;
1073:       bjac->starts[i] = start;

1075:       ISCreateStride(PETSC_COMM_SELF,m,start,1,&is);
1076:       bjac->is[i] = is;
1077:       PetscLogObjectParent((PetscObject)pc,(PetscObject)is);

1079:       start += m;
1080:     }
1081:   } else {
1082:     bjac = (PC_BJacobi_Multiblock*)jac->data;
1083:     /*
1084:        Destroy the blocks from the previous iteration
1085:     */
1086:     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
1087:       MatDestroyMatrices(n_local,&bjac->pmat);
1088:       if (pc->useAmat) {
1089:         MatDestroyMatrices(n_local,&bjac->mat);
1090:       }
1091:       scall = MAT_INITIAL_MATRIX;
1092:     } else scall = MAT_REUSE_MATRIX;
1093:   }

1095:   MatCreateSubMatrices(pmat,n_local,bjac->is,bjac->is,scall,&bjac->pmat);
1096:   if (pc->useAmat) {
1097:     MatCreateSubMatrices(mat,n_local,bjac->is,bjac->is,scall,&bjac->mat);
1098:   }
1099:   /* Return control to the user so that the submatrices can be modified (e.g., to apply
1100:      different boundary conditions for the submatrices than for the global problem) */
1101:   PCModifySubMatrices(pc,n_local,bjac->is,bjac->is,bjac->pmat,pc->modifysubmatricesP);

1103:   PetscObjectGetOptionsPrefix((PetscObject)pmat,&pprefix);
1104:   for (i=0; i<n_local; i++) {
1105:     PetscLogObjectParent((PetscObject)pc,(PetscObject)bjac->pmat[i]);
1106:     PetscObjectSetOptionsPrefix((PetscObject)bjac->pmat[i],pprefix);
1107:     if (pc->useAmat) {
1108:       PetscLogObjectParent((PetscObject)pc,(PetscObject)bjac->mat[i]);
1109:       PetscObjectGetOptionsPrefix((PetscObject)mat,&mprefix);
1110:       PetscObjectSetOptionsPrefix((PetscObject)bjac->mat[i],mprefix);
1111:       KSPSetOperators(jac->ksp[i],bjac->mat[i],bjac->pmat[i]);
1112:     } else {
1113:       KSPSetOperators(jac->ksp[i],bjac->pmat[i],bjac->pmat[i]);
1114:     }
1115:     if (pc->setfromoptionscalled) {
1116:       KSPSetFromOptions(jac->ksp[i]);
1117:     }
1118:   }
1119:   return(0);
1120: }

1122: /* ---------------------------------------------------------------------------------------------*/
1123: /*
1124:       These are for a single block with multiple processes;
1125: */
1126: static PetscErrorCode PCReset_BJacobi_Multiproc(PC pc)
1127: {
1128:   PC_BJacobi           *jac   = (PC_BJacobi*)pc->data;
1129:   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1130:   PetscErrorCode       ierr;

1133:   VecDestroy(&mpjac->ysub);
1134:   VecDestroy(&mpjac->xsub);
1135:   MatDestroy(&mpjac->submats);
1136:   if (jac->ksp) {KSPReset(jac->ksp[0]);}
1137:   return(0);
1138: }

1140: static PetscErrorCode PCDestroy_BJacobi_Multiproc(PC pc)
1141: {
1142:   PC_BJacobi           *jac   = (PC_BJacobi*)pc->data;
1143:   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1144:   PetscErrorCode       ierr;

1147:   PCReset_BJacobi_Multiproc(pc);
1148:   KSPDestroy(&jac->ksp[0]);
1149:   PetscFree(jac->ksp);
1150:   PetscSubcommDestroy(&mpjac->psubcomm);

1152:   PetscFree(mpjac);
1153:   PetscFree(pc->data);
1154:   return(0);
1155: }

1157: static PetscErrorCode PCApply_BJacobi_Multiproc(PC pc,Vec x,Vec y)
1158: {
1159:   PC_BJacobi           *jac   = (PC_BJacobi*)pc->data;
1160:   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1161:   PetscErrorCode       ierr;
1162:   PetscScalar          *yarray;
1163:   const PetscScalar    *xarray;
1164:   KSPConvergedReason   reason;

1167:   /* place x's and y's local arrays into xsub and ysub */
1168:   VecGetArrayRead(x,&xarray);
1169:   VecGetArray(y,&yarray);
1170:   VecPlaceArray(mpjac->xsub,xarray);
1171:   VecPlaceArray(mpjac->ysub,yarray);

1173:   /* apply preconditioner on each matrix block */
1174:   PetscLogEventBegin(PC_ApplyOnBlocks,jac->ksp[0],mpjac->xsub,mpjac->ysub,0);
1175:   KSPSolve(jac->ksp[0],mpjac->xsub,mpjac->ysub);
1176:   KSPCheckSolve(jac->ksp[0],pc,mpjac->ysub);
1177:   PetscLogEventEnd(PC_ApplyOnBlocks,jac->ksp[0],mpjac->xsub,mpjac->ysub,0);
1178:   KSPGetConvergedReason(jac->ksp[0],&reason);
1179:   if (reason == KSP_DIVERGED_PC_FAILED) {
1180:     pc->failedreason = PC_SUBPC_ERROR;
1181:   }

1183:   VecResetArray(mpjac->xsub);
1184:   VecResetArray(mpjac->ysub);
1185:   VecRestoreArrayRead(x,&xarray);
1186:   VecRestoreArray(y,&yarray);
1187:   return(0);
1188: }

1190: static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC pc)
1191: {
1192:   PC_BJacobi           *jac   = (PC_BJacobi*)pc->data;
1193:   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1194:   PetscErrorCode       ierr;
1195:   PetscInt             m,n;
1196:   MPI_Comm             comm,subcomm=0;
1197:   const char           *prefix;
1198:   PetscBool            wasSetup = PETSC_TRUE;
1199: #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_VIENNACL)
1200:   PetscBool              is_gpumatrix = PETSC_FALSE;
1201: #endif


1205:   PetscObjectGetComm((PetscObject)pc,&comm);
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:     PetscNewLog(pc,&mpjac);
1211:     jac->data = (void*)mpjac;

1213:     /* initialize datastructure mpjac */
1214:     if (!jac->psubcomm) {
1215:       /* Create default contiguous subcommunicatiors if user does not provide them */
1216:       PetscSubcommCreate(comm,&jac->psubcomm);
1217:       PetscSubcommSetNumber(jac->psubcomm,jac->n);
1218:       PetscSubcommSetType(jac->psubcomm,PETSC_SUBCOMM_CONTIGUOUS);
1219:       PetscLogObjectMemory((PetscObject)pc,sizeof(PetscSubcomm));
1220:     }
1221:     mpjac->psubcomm = jac->psubcomm;
1222:     subcomm         = PetscSubcommChild(mpjac->psubcomm);

1224:     /* Get matrix blocks of pmat */
1225:     MatGetMultiProcBlock(pc->pmat,subcomm,MAT_INITIAL_MATRIX,&mpjac->submats);

1227:     /* create a new PC that processors in each subcomm have copy of */
1228:     PetscMalloc1(1,&jac->ksp);
1229:     KSPCreate(subcomm,&jac->ksp[0]);
1230:     KSPSetErrorIfNotConverged(jac->ksp[0],pc->erroriffailure);
1231:     PetscObjectIncrementTabLevel((PetscObject)jac->ksp[0],(PetscObject)pc,1);
1232:     PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->ksp[0]);
1233:     KSPSetOperators(jac->ksp[0],mpjac->submats,mpjac->submats);
1234:     KSPGetPC(jac->ksp[0],&mpjac->pc);

1236:     PCGetOptionsPrefix(pc,&prefix);
1237:     KSPSetOptionsPrefix(jac->ksp[0],prefix);
1238:     KSPAppendOptionsPrefix(jac->ksp[0],"sub_");
1239:     /*
1240:       PetscMPIInt rank,subsize,subrank;
1241:       MPI_Comm_rank(comm,&rank);
1242:       MPI_Comm_size(subcomm,&subsize);
1243:       MPI_Comm_rank(subcomm,&subrank);

1245:       MatGetLocalSize(mpjac->submats,&m,NULL);
1246:       MatGetSize(mpjac->submats,&n,NULL);
1247:       PetscSynchronizedPrintf(comm,"[%d], sub-size %d,sub-rank %d\n",rank,subsize,subrank);
1248:       PetscSynchronizedFlush(comm,PETSC_STDOUT);
1249:     */

1251:     /* create dummy vectors xsub and ysub */
1252:     MatGetLocalSize(mpjac->submats,&m,&n);
1253:     VecCreateMPIWithArray(subcomm,1,n,PETSC_DECIDE,NULL,&mpjac->xsub);
1254:     VecCreateMPIWithArray(subcomm,1,m,PETSC_DECIDE,NULL,&mpjac->ysub);
1255: #if defined(PETSC_HAVE_CUDA)
1256:     PetscObjectTypeCompareAny((PetscObject)mpjac->submats,&is_gpumatrix,MATAIJCUSPARSE,MATMPIAIJCUSPARSE,"");
1257:     if (is_gpumatrix) {
1258:       VecSetType(mpjac->xsub,VECMPICUDA);
1259:       VecSetType(mpjac->ysub,VECMPICUDA);
1260:     }
1261: #endif
1262: #if defined(PETSC_HAVE_VIENNACL)
1263:     PetscObjectTypeCompareAny((PetscObject)mpjac->submats,&is_gpumatrix,MATAIJVIENNACL,MATMPIAIJVIENNACL,"");
1264:     if (is_gpumatrix) {
1265:       VecSetType(mpjac->xsub,VECMPIVIENNACL);
1266:       VecSetType(mpjac->ysub,VECMPIVIENNACL);
1267:     }
1268: #endif
1269:     PetscLogObjectParent((PetscObject)pc,(PetscObject)mpjac->xsub);
1270:     PetscLogObjectParent((PetscObject)pc,(PetscObject)mpjac->ysub);

1272:     pc->ops->reset   = PCReset_BJacobi_Multiproc;
1273:     pc->ops->destroy = PCDestroy_BJacobi_Multiproc;
1274:     pc->ops->apply   = PCApply_BJacobi_Multiproc;
1275:   } else { /* pc->setupcalled */
1276:     subcomm = PetscSubcommChild(mpjac->psubcomm);
1277:     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
1278:       /* destroy old matrix blocks, then get new matrix blocks */
1279:       if (mpjac->submats) {MatDestroy(&mpjac->submats);}
1280:       MatGetMultiProcBlock(pc->pmat,subcomm,MAT_INITIAL_MATRIX,&mpjac->submats);
1281:     } else {
1282:       MatGetMultiProcBlock(pc->pmat,subcomm,MAT_REUSE_MATRIX,&mpjac->submats);
1283:     }
1284:     KSPSetOperators(jac->ksp[0],mpjac->submats,mpjac->submats);
1285:   }

1287:   if (!wasSetup && pc->setfromoptionscalled) {
1288:     KSPSetFromOptions(jac->ksp[0]);
1289:   }
1290:   return(0);
1291: }