Actual source code: bjacobi.c

petsc-3.12.5 2020-03-29
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:         PetscArraycpy(jac->l_lens,jac->g_lens,jac->n_local);
 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:         PetscArraycpy(jac->l_lens,jac->g_lens+i_start,jac->n_local);
 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:     PetscArraycpy(jac->g_lens,lens,blocks);
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:     PetscArraycpy(jac->l_lens,lens,blocks);
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: .seealso: PCBJacobiGetSubKSP()
367: @*/
368: PetscErrorCode  PCBJacobiGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
369: {

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

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

382:    Collective on PC

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

389:    Options Database Key:
390: .  -pc_bjacobi_blocks <blocks> - Sets the number of global blocks

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

396:    Level: intermediate

398: .seealso: PCSetUseAmat(), PCBJacobiSetLocalBlocks()
399: @*/
400: PetscErrorCode  PCBJacobiSetTotalBlocks(PC pc,PetscInt blocks,const PetscInt lens[])
401: {

406:   if (blocks <= 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Must have positive blocks");
407:   PetscTryMethod(pc,"PCBJacobiSetTotalBlocks_C",(PC,PetscInt,const PetscInt[]),(pc,blocks,lens));
408:   return(0);
409: }

411: /*@C
412:    PCBJacobiGetTotalBlocks - Gets the global number of blocks for the block
413:    Jacobi preconditioner.

415:    Not Collective

417:    Input Parameter:
418: .  pc - the preconditioner context

420:    Output parameters:
421: +  blocks - the number of blocks
422: -  lens - integer array containing the size of each block

424:    Level: intermediate

426: .seealso: PCSetUseAmat(), PCBJacobiGetLocalBlocks()
427: @*/
428: PetscErrorCode  PCBJacobiGetTotalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
429: {

435:   PetscUseMethod(pc,"PCBJacobiGetTotalBlocks_C",(PC,PetscInt*, const PetscInt *[]),(pc,blocks,lens));
436:   return(0);
437: }

439: /*@
440:    PCBJacobiSetLocalBlocks - Sets the local number of blocks for the block
441:    Jacobi preconditioner.

443:    Not Collective

445:    Input Parameters:
446: +  pc - the preconditioner context
447: .  blocks - the number of blocks
448: -  lens - [optional] integer array containing size of each block

450:    Options Database Key:
451: .  -pc_bjacobi_local_blocks <blocks> - Sets the number of local blocks

453:    Note:
454:    Currently only a limited number of blocking configurations are supported.

456:    Level: intermediate

458: .seealso: PCSetUseAmat(), PCBJacobiSetTotalBlocks()
459: @*/
460: PetscErrorCode  PCBJacobiSetLocalBlocks(PC pc,PetscInt blocks,const PetscInt lens[])
461: {

466:   if (blocks < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have nonegative blocks");
467:   PetscTryMethod(pc,"PCBJacobiSetLocalBlocks_C",(PC,PetscInt,const PetscInt []),(pc,blocks,lens));
468:   return(0);
469: }

471: /*@C
472:    PCBJacobiGetLocalBlocks - Gets the local number of blocks for the block
473:    Jacobi preconditioner.

475:    Not Collective

477:    Input Parameters:
478: +  pc - the preconditioner context
479: .  blocks - the number of blocks
480: -  lens - [optional] integer array containing size of each block

482:    Note:
483:    Currently only a limited number of blocking configurations are supported.

485:    Level: intermediate

487: .seealso: PCSetUseAmat(), PCBJacobiGetTotalBlocks()
488: @*/
489: PetscErrorCode  PCBJacobiGetLocalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
490: {

496:   PetscUseMethod(pc,"PCBJacobiGetLocalBlocks_C",(PC,PetscInt*, const PetscInt *[]),(pc,blocks,lens));
497:   return(0);
498: }

500: /* -----------------------------------------------------------------------------------*/

502: /*MC
503:    PCBJACOBI - Use block Jacobi preconditioning, each block is (approximately) solved with
504:            its own KSP object.

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

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

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

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

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

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

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

528:    Level: beginner

530: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
531:            PCASM, PCSetUseAmat(), PCGetUseAmat(), PCBJacobiGetSubKSP(), PCBJacobiSetTotalBlocks(),
532:            PCBJacobiSetLocalBlocks(), PCSetModifySubMatrices()
533: M*/

535: PETSC_EXTERN PetscErrorCode PCCreate_BJacobi(PC pc)
536: {
538:   PetscMPIInt    rank;
539:   PC_BJacobi     *jac;

542:   PetscNewLog(pc,&jac);
543:   MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);

545:   pc->ops->apply           = 0;
546:   pc->ops->applytranspose  = 0;
547:   pc->ops->setup           = PCSetUp_BJacobi;
548:   pc->ops->destroy         = PCDestroy_BJacobi;
549:   pc->ops->setfromoptions  = PCSetFromOptions_BJacobi;
550:   pc->ops->view            = PCView_BJacobi;
551:   pc->ops->applyrichardson = 0;

553:   pc->data               = (void*)jac;
554:   jac->n                 = -1;
555:   jac->n_local           = -1;
556:   jac->first_local       = rank;
557:   jac->ksp               = 0;
558:   jac->same_local_solves = PETSC_TRUE;
559:   jac->g_lens            = 0;
560:   jac->l_lens            = 0;
561:   jac->psubcomm          = 0;

563:   PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiGetSubKSP_C",PCBJacobiGetSubKSP_BJacobi);
564:   PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiSetTotalBlocks_C",PCBJacobiSetTotalBlocks_BJacobi);
565:   PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiGetTotalBlocks_C",PCBJacobiGetTotalBlocks_BJacobi);
566:   PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiSetLocalBlocks_C",PCBJacobiSetLocalBlocks_BJacobi);
567:   PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiGetLocalBlocks_C",PCBJacobiGetLocalBlocks_BJacobi);
568:   return(0);
569: }

571: /* --------------------------------------------------------------------------------------------*/
572: /*
573:         These are for a single block per processor; works for AIJ, BAIJ; Seq and MPI
574: */
575: static PetscErrorCode PCReset_BJacobi_Singleblock(PC pc)
576: {
577:   PC_BJacobi             *jac  = (PC_BJacobi*)pc->data;
578:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
579:   PetscErrorCode         ierr;

582:   KSPReset(jac->ksp[0]);
583:   VecDestroy(&bjac->x);
584:   VecDestroy(&bjac->y);
585:   return(0);
586: }

588: static PetscErrorCode PCDestroy_BJacobi_Singleblock(PC pc)
589: {
590:   PC_BJacobi             *jac  = (PC_BJacobi*)pc->data;
591:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
592:   PetscErrorCode         ierr;

595:   PCReset_BJacobi_Singleblock(pc);
596:   KSPDestroy(&jac->ksp[0]);
597:   PetscFree(jac->ksp);
598:   PetscFree(jac->l_lens);
599:   PetscFree(jac->g_lens);
600:   PetscFree(bjac);
601:   PetscFree(pc->data);
602:   return(0);
603: }

605: static PetscErrorCode PCSetUpOnBlocks_BJacobi_Singleblock(PC pc)
606: {
607:   PetscErrorCode     ierr;
608:   PC_BJacobi         *jac = (PC_BJacobi*)pc->data;
609:   KSP                subksp = jac->ksp[0];
610:   KSPConvergedReason reason;

613:   KSPSetUp(subksp);
614:   KSPGetConvergedReason(subksp,&reason);
615:   if (reason == KSP_DIVERGED_PC_FAILED) {
616:     pc->failedreason = PC_SUBPC_ERROR;
617:   }
618:   return(0);
619: }

621: static PetscErrorCode PCApply_BJacobi_Singleblock(PC pc,Vec x,Vec y)
622: {
623:   PetscErrorCode         ierr;
624:   PC_BJacobi             *jac  = (PC_BJacobi*)pc->data;
625:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;

628:   VecGetLocalVectorRead(x, bjac->x);
629:   VecGetLocalVector(y, bjac->y);
630:  /* Since the inner KSP matrix may point directly to the diagonal block of an MPI matrix the inner
631:      matrix may change even if the outter KSP/PC has not updated the preconditioner, this will trigger a rebuild
632:      of the inner preconditioner automatically unless we pass down the outter preconditioners reuse flag.*/
633:   KSPSetReusePreconditioner(jac->ksp[0],pc->reusepreconditioner);
634:   KSPSolve(jac->ksp[0],bjac->x,bjac->y);
635:   KSPCheckSolve(jac->ksp[0],pc,bjac->y);
636:   VecRestoreLocalVectorRead(x, bjac->x);
637:   VecRestoreLocalVector(y, bjac->y);
638:   return(0);
639: }

641: static PetscErrorCode PCApplySymmetricLeft_BJacobi_Singleblock(PC pc,Vec x,Vec y)
642: {
643:   PetscErrorCode         ierr;
644:   PC_BJacobi             *jac  = (PC_BJacobi*)pc->data;
645:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
646:   PetscScalar            *y_array;
647:   const PetscScalar      *x_array;
648:   PC                     subpc;

651:   /*
652:       The VecPlaceArray() is to avoid having to copy the
653:     y vector into the bjac->x vector. The reason for
654:     the bjac->x vector is that we need a sequential vector
655:     for the sequential solve.
656:   */
657:   VecGetArrayRead(x,&x_array);
658:   VecGetArray(y,&y_array);
659:   VecPlaceArray(bjac->x,x_array);
660:   VecPlaceArray(bjac->y,y_array);
661:   /* apply the symmetric left portion of the inner PC operator */
662:   /* note this by-passes the inner KSP and its options completely */
663:   KSPGetPC(jac->ksp[0],&subpc);
664:   PCApplySymmetricLeft(subpc,bjac->x,bjac->y);
665:   VecResetArray(bjac->x);
666:   VecResetArray(bjac->y);
667:   VecRestoreArrayRead(x,&x_array);
668:   VecRestoreArray(y,&y_array);
669:   return(0);
670: }

672: static PetscErrorCode PCApplySymmetricRight_BJacobi_Singleblock(PC pc,Vec x,Vec y)
673: {
674:   PetscErrorCode         ierr;
675:   PC_BJacobi             *jac  = (PC_BJacobi*)pc->data;
676:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
677:   PetscScalar            *y_array;
678:   const PetscScalar      *x_array;
679:   PC                     subpc;

682:   /*
683:       The VecPlaceArray() is to avoid having to copy the
684:     y vector into the bjac->x vector. The reason for
685:     the bjac->x vector is that we need a sequential vector
686:     for the sequential solve.
687:   */
688:   VecGetArrayRead(x,&x_array);
689:   VecGetArray(y,&y_array);
690:   VecPlaceArray(bjac->x,x_array);
691:   VecPlaceArray(bjac->y,y_array);

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

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

699:   VecRestoreArrayRead(x,&x_array);
700:   VecRestoreArray(y,&y_array);
701:   return(0);
702: }

704: static PetscErrorCode PCApplyTranspose_BJacobi_Singleblock(PC pc,Vec x,Vec y)
705: {
706:   PetscErrorCode         ierr;
707:   PC_BJacobi             *jac  = (PC_BJacobi*)pc->data;
708:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
709:   PetscScalar            *y_array;
710:   const PetscScalar      *x_array;

713:   /*
714:       The VecPlaceArray() is to avoid having to copy the
715:     y vector into the bjac->x vector. The reason for
716:     the bjac->x vector is that we need a sequential vector
717:     for the sequential solve.
718:   */
719:   VecGetArrayRead(x,&x_array);
720:   VecGetArray(y,&y_array);
721:   VecPlaceArray(bjac->x,x_array);
722:   VecPlaceArray(bjac->y,y_array);
723:   KSPSolveTranspose(jac->ksp[0],bjac->x,bjac->y);
724:   KSPCheckSolve(jac->ksp[0],pc,bjac->y);
725:   VecResetArray(bjac->x);
726:   VecResetArray(bjac->y);
727:   VecRestoreArrayRead(x,&x_array);
728:   VecRestoreArray(y,&y_array);
729:   return(0);
730: }

732: static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC pc,Mat mat,Mat pmat)
733: {
734:   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
735:   PetscErrorCode         ierr;
736:   PetscInt               m;
737:   KSP                    ksp;
738:   PC_BJacobi_Singleblock *bjac;
739:   PetscBool              wasSetup = PETSC_TRUE;
740: #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_VIENNACL)
741:   PetscBool              is_gpumatrix = PETSC_FALSE;
742: #endif

745:   if (!pc->setupcalled) {
746:     const char *prefix;

748:     if (!jac->ksp) {
749:       wasSetup = PETSC_FALSE;

751:       KSPCreate(PETSC_COMM_SELF,&ksp);
752:       KSPSetErrorIfNotConverged(ksp,pc->erroriffailure);
753:       PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);
754:       PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);
755:       KSPSetType(ksp,KSPPREONLY);
756:       PCGetOptionsPrefix(pc,&prefix);
757:       KSPSetOptionsPrefix(ksp,prefix);
758:       KSPAppendOptionsPrefix(ksp,"sub_");

760:       pc->ops->reset               = PCReset_BJacobi_Singleblock;
761:       pc->ops->destroy             = PCDestroy_BJacobi_Singleblock;
762:       pc->ops->apply               = PCApply_BJacobi_Singleblock;
763:       pc->ops->applysymmetricleft  = PCApplySymmetricLeft_BJacobi_Singleblock;
764:       pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Singleblock;
765:       pc->ops->applytranspose      = PCApplyTranspose_BJacobi_Singleblock;
766:       pc->ops->setuponblocks       = PCSetUpOnBlocks_BJacobi_Singleblock;

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

771:       PetscNewLog(pc,&bjac);
772:       jac->data = (void*)bjac;
773:     } else {
774:       ksp  = jac->ksp[0];
775:       bjac = (PC_BJacobi_Singleblock*)jac->data;
776:     }

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

820: /* ---------------------------------------------------------------------------------------------*/
821: static PetscErrorCode PCReset_BJacobi_Multiblock(PC pc)
822: {
823:   PC_BJacobi            *jac  = (PC_BJacobi*)pc->data;
824:   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
825:   PetscErrorCode        ierr;
826:   PetscInt              i;

829:   if (bjac && bjac->pmat) {
830:     MatDestroyMatrices(jac->n_local,&bjac->pmat);
831:     if (pc->useAmat) {
832:       MatDestroyMatrices(jac->n_local,&bjac->mat);
833:     }
834:   }

836:   for (i=0; i<jac->n_local; i++) {
837:     KSPReset(jac->ksp[i]);
838:     if (bjac && bjac->x) {
839:       VecDestroy(&bjac->x[i]);
840:       VecDestroy(&bjac->y[i]);
841:       ISDestroy(&bjac->is[i]);
842:     }
843:   }
844:   PetscFree(jac->l_lens);
845:   PetscFree(jac->g_lens);
846:   return(0);
847: }

849: static PetscErrorCode PCDestroy_BJacobi_Multiblock(PC pc)
850: {
851:   PC_BJacobi            *jac  = (PC_BJacobi*)pc->data;
852:   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
853:   PetscErrorCode        ierr;
854:   PetscInt              i;

857:   PCReset_BJacobi_Multiblock(pc);
858:   if (bjac) {
859:     PetscFree2(bjac->x,bjac->y);
860:     PetscFree(bjac->starts);
861:     PetscFree(bjac->is);
862:   }
863:   PetscFree(jac->data);
864:   for (i=0; i<jac->n_local; i++) {
865:     KSPDestroy(&jac->ksp[i]);
866:   }
867:   PetscFree(jac->ksp);
868:   PetscFree(pc->data);
869:   return(0);
870: }

872: static PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiblock(PC pc)
873: {
874:   PC_BJacobi         *jac = (PC_BJacobi*)pc->data;
875:   PetscErrorCode     ierr;
876:   PetscInt           i,n_local = jac->n_local;
877:   KSPConvergedReason reason;

880:   for (i=0; i<n_local; i++) {
881:     KSPSetUp(jac->ksp[i]);
882:     KSPGetConvergedReason(jac->ksp[i],&reason);
883:     if (reason == KSP_DIVERGED_PC_FAILED) {
884:       pc->failedreason = PC_SUBPC_ERROR;
885:     }
886:   }
887:   return(0);
888: }

890: /*
891:       Preconditioner for block Jacobi
892: */
893: static PetscErrorCode PCApply_BJacobi_Multiblock(PC pc,Vec x,Vec y)
894: {
895:   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
896:   PetscErrorCode        ierr;
897:   PetscInt              i,n_local = jac->n_local;
898:   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
899:   PetscScalar           *yin;
900:   const PetscScalar     *xin;

903:   VecGetArrayRead(x,&xin);
904:   VecGetArray(y,&yin);
905:   for (i=0; i<n_local; i++) {
906:     /*
907:        To avoid copying the subvector from x into a workspace we instead
908:        make the workspace vector array point to the subpart of the array of
909:        the global vector.
910:     */
911:     VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);
912:     VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);

914:     PetscLogEventBegin(PC_ApplyOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);
915:     KSPSolve(jac->ksp[i],bjac->x[i],bjac->y[i]);
916:     KSPCheckSolve(jac->ksp[i],pc,bjac->y[i]);
917:     PetscLogEventEnd(PC_ApplyOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);

919:     VecResetArray(bjac->x[i]);
920:     VecResetArray(bjac->y[i]);
921:   }
922:   VecRestoreArrayRead(x,&xin);
923:   VecRestoreArray(y,&yin);
924:   return(0);
925: }

927: /*
928:       Preconditioner for block Jacobi
929: */
930: static PetscErrorCode PCApplyTranspose_BJacobi_Multiblock(PC pc,Vec x,Vec y)
931: {
932:   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
933:   PetscErrorCode        ierr;
934:   PetscInt              i,n_local = jac->n_local;
935:   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
936:   PetscScalar           *yin;
937:   const PetscScalar     *xin;

940:   VecGetArrayRead(x,&xin);
941:   VecGetArray(y,&yin);
942:   for (i=0; i<n_local; i++) {
943:     /*
944:        To avoid copying the subvector from x into a workspace we instead
945:        make the workspace vector array point to the subpart of the array of
946:        the global vector.
947:     */
948:     VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);
949:     VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);

951:     PetscLogEventBegin(PC_ApplyTransposeOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);
952:     KSPSolveTranspose(jac->ksp[i],bjac->x[i],bjac->y[i]);
953:     KSPCheckSolve(jac->ksp[i],pc,bjac->y[i]);
954:     PetscLogEventEnd(PC_ApplyTransposeOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);

956:     VecResetArray(bjac->x[i]);
957:     VecResetArray(bjac->y[i]);
958:   }
959:   VecRestoreArrayRead(x,&xin);
960:   VecRestoreArray(y,&yin);
961:   return(0);
962: }

964: static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC pc,Mat mat,Mat pmat)
965: {
966:   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
967:   PetscErrorCode        ierr;
968:   PetscInt              m,n_local,N,M,start,i;
969:   const char            *prefix,*pprefix,*mprefix;
970:   KSP                   ksp;
971:   Vec                   x,y;
972:   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
973:   PC                    subpc;
974:   IS                    is;
975:   MatReuse              scall;
976: #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_VIENNACL)
977:   PetscBool              is_gpumatrix = PETSC_FALSE;
978: #endif

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

983:   n_local = jac->n_local;

985:   if (pc->useAmat) {
986:     PetscBool same;
987:     PetscObjectTypeCompare((PetscObject)mat,((PetscObject)pmat)->type_name,&same);
988:     if (!same) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Matrices not of same type");
989:   }

991:   if (!pc->setupcalled) {
992:     scall = MAT_INITIAL_MATRIX;

994:     if (!jac->ksp) {
995:       pc->ops->reset         = PCReset_BJacobi_Multiblock;
996:       pc->ops->destroy       = PCDestroy_BJacobi_Multiblock;
997:       pc->ops->apply         = PCApply_BJacobi_Multiblock;
998:       pc->ops->applytranspose= PCApplyTranspose_BJacobi_Multiblock;
999:       pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiblock;

1001:       PetscNewLog(pc,&bjac);
1002:       PetscMalloc1(n_local,&jac->ksp);
1003:       PetscLogObjectMemory((PetscObject)pc,sizeof(n_local*sizeof(KSP)));
1004:       PetscMalloc2(n_local,&bjac->x,n_local,&bjac->y);
1005:       PetscMalloc1(n_local,&bjac->starts);
1006:       PetscLogObjectMemory((PetscObject)pc,sizeof(n_local*sizeof(PetscScalar)));

1008:       jac->data = (void*)bjac;
1009:       PetscMalloc1(n_local,&bjac->is);
1010:       PetscLogObjectMemory((PetscObject)pc,sizeof(n_local*sizeof(IS)));

1012:       for (i=0; i<n_local; i++) {
1013:         KSPCreate(PETSC_COMM_SELF,&ksp);
1014:         KSPSetErrorIfNotConverged(ksp,pc->erroriffailure);
1015:         PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);
1016:         PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);
1017:         KSPSetType(ksp,KSPPREONLY);
1018:         KSPGetPC(ksp,&subpc);
1019:         PCGetOptionsPrefix(pc,&prefix);
1020:         KSPSetOptionsPrefix(ksp,prefix);
1021:         KSPAppendOptionsPrefix(ksp,"sub_");

1023:         jac->ksp[i] = ksp;
1024:       }
1025:     } else {
1026:       bjac = (PC_BJacobi_Multiblock*)jac->data;
1027:     }

1029:     start = 0;
1030:     for (i=0; i<n_local; i++) {
1031:       m = jac->l_lens[i];
1032:       /*
1033:       The reason we need to generate these vectors is to serve
1034:       as the right-hand side and solution vector for the solve on the
1035:       block. We do not need to allocate space for the vectors since
1036:       that is provided via VecPlaceArray() just before the call to
1037:       KSPSolve() on the block.

1039:       */
1040:       VecCreateSeq(PETSC_COMM_SELF,m,&x);
1041:       VecCreateSeqWithArray(PETSC_COMM_SELF,1,m,NULL,&y);
1042: #if defined(PETSC_HAVE_CUDA)
1043:       PetscObjectTypeCompareAny((PetscObject)pmat,&is_gpumatrix,MATAIJCUSPARSE,MATSEQAIJCUSPARSE,MATMPIAIJCUSPARSE,"");
1044:       if (is_gpumatrix) {
1045:         VecSetType(x,VECCUDA);
1046:         VecSetType(y,VECCUDA);
1047:       }
1048: #endif
1049: #if defined(PETSC_HAVE_VIENNACL)
1050:       PetscObjectTypeCompareAny((PetscObject)pmat,&is_gpumatrix,MATAIJVIENNACL,MATSEQAIJVIENNACL,MATMPIAIJVIENNACL,"");
1051:       if (is_gpumatrix) {
1052:         VecSetType(x,VECVIENNACL);
1053:         VecSetType(y,VECVIENNACL);
1054:       }
1055: #endif
1056:       PetscLogObjectParent((PetscObject)pc,(PetscObject)x);
1057:       PetscLogObjectParent((PetscObject)pc,(PetscObject)y);

1059:       bjac->x[i]      = x;
1060:       bjac->y[i]      = y;
1061:       bjac->starts[i] = start;

1063:       ISCreateStride(PETSC_COMM_SELF,m,start,1,&is);
1064:       bjac->is[i] = is;
1065:       PetscLogObjectParent((PetscObject)pc,(PetscObject)is);

1067:       start += m;
1068:     }
1069:   } else {
1070:     bjac = (PC_BJacobi_Multiblock*)jac->data;
1071:     /*
1072:        Destroy the blocks from the previous iteration
1073:     */
1074:     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
1075:       MatDestroyMatrices(n_local,&bjac->pmat);
1076:       if (pc->useAmat) {
1077:         MatDestroyMatrices(n_local,&bjac->mat);
1078:       }
1079:       scall = MAT_INITIAL_MATRIX;
1080:     } else scall = MAT_REUSE_MATRIX;
1081:   }

1083:   MatCreateSubMatrices(pmat,n_local,bjac->is,bjac->is,scall,&bjac->pmat);
1084:   if (pc->useAmat) {
1085:     MatCreateSubMatrices(mat,n_local,bjac->is,bjac->is,scall,&bjac->mat);
1086:   }
1087:   /* Return control to the user so that the submatrices can be modified (e.g., to apply
1088:      different boundary conditions for the submatrices than for the global problem) */
1089:   PCModifySubMatrices(pc,n_local,bjac->is,bjac->is,bjac->pmat,pc->modifysubmatricesP);

1091:   PetscObjectGetOptionsPrefix((PetscObject)pmat,&pprefix);
1092:   for (i=0; i<n_local; i++) {
1093:     PetscLogObjectParent((PetscObject)pc,(PetscObject)bjac->pmat[i]);
1094:     PetscObjectSetOptionsPrefix((PetscObject)bjac->pmat[i],pprefix);
1095:     if (pc->useAmat) {
1096:       PetscLogObjectParent((PetscObject)pc,(PetscObject)bjac->mat[i]);
1097:       PetscObjectGetOptionsPrefix((PetscObject)mat,&mprefix);
1098:       PetscObjectSetOptionsPrefix((PetscObject)bjac->mat[i],mprefix);
1099:       KSPSetOperators(jac->ksp[i],bjac->mat[i],bjac->pmat[i]);
1100:     } else {
1101:       KSPSetOperators(jac->ksp[i],bjac->pmat[i],bjac->pmat[i]);
1102:     }
1103:     if (pc->setfromoptionscalled) {
1104:       KSPSetFromOptions(jac->ksp[i]);
1105:     }
1106:   }
1107:   return(0);
1108: }

1110: /* ---------------------------------------------------------------------------------------------*/
1111: /*
1112:       These are for a single block with multiple processes;
1113: */
1114: static PetscErrorCode PCReset_BJacobi_Multiproc(PC pc)
1115: {
1116:   PC_BJacobi           *jac   = (PC_BJacobi*)pc->data;
1117:   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1118:   PetscErrorCode       ierr;

1121:   VecDestroy(&mpjac->ysub);
1122:   VecDestroy(&mpjac->xsub);
1123:   MatDestroy(&mpjac->submats);
1124:   if (jac->ksp) {KSPReset(jac->ksp[0]);}
1125:   return(0);
1126: }

1128: static PetscErrorCode PCDestroy_BJacobi_Multiproc(PC pc)
1129: {
1130:   PC_BJacobi           *jac   = (PC_BJacobi*)pc->data;
1131:   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1132:   PetscErrorCode       ierr;

1135:   PCReset_BJacobi_Multiproc(pc);
1136:   KSPDestroy(&jac->ksp[0]);
1137:   PetscFree(jac->ksp);
1138:   PetscSubcommDestroy(&mpjac->psubcomm);

1140:   PetscFree(mpjac);
1141:   PetscFree(pc->data);
1142:   return(0);
1143: }

1145: static PetscErrorCode PCApply_BJacobi_Multiproc(PC pc,Vec x,Vec y)
1146: {
1147:   PC_BJacobi           *jac   = (PC_BJacobi*)pc->data;
1148:   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1149:   PetscErrorCode       ierr;
1150:   PetscScalar          *yarray;
1151:   const PetscScalar    *xarray;
1152:   KSPConvergedReason   reason;

1155:   /* place x's and y's local arrays into xsub and ysub */
1156:   VecGetArrayRead(x,&xarray);
1157:   VecGetArray(y,&yarray);
1158:   VecPlaceArray(mpjac->xsub,xarray);
1159:   VecPlaceArray(mpjac->ysub,yarray);

1161:   /* apply preconditioner on each matrix block */
1162:   PetscLogEventBegin(PC_ApplyOnBlocks,jac->ksp[0],mpjac->xsub,mpjac->ysub,0);
1163:   KSPSolve(jac->ksp[0],mpjac->xsub,mpjac->ysub);
1164:   KSPCheckSolve(jac->ksp[0],pc,mpjac->ysub);
1165:   PetscLogEventEnd(PC_ApplyOnBlocks,jac->ksp[0],mpjac->xsub,mpjac->ysub,0);
1166:   KSPGetConvergedReason(jac->ksp[0],&reason);
1167:   if (reason == KSP_DIVERGED_PC_FAILED) {
1168:     pc->failedreason = PC_SUBPC_ERROR;
1169:   }

1171:   VecResetArray(mpjac->xsub);
1172:   VecResetArray(mpjac->ysub);
1173:   VecRestoreArrayRead(x,&xarray);
1174:   VecRestoreArray(y,&yarray);
1175:   return(0);
1176: }

1178: static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC pc)
1179: {
1180:   PC_BJacobi           *jac   = (PC_BJacobi*)pc->data;
1181:   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1182:   PetscErrorCode       ierr;
1183:   PetscInt             m,n;
1184:   MPI_Comm             comm,subcomm=0;
1185:   const char           *prefix;
1186:   PetscBool            wasSetup = PETSC_TRUE;
1187: #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_VIENNACL)
1188:   PetscBool              is_gpumatrix = PETSC_FALSE;
1189: #endif


1193:   PetscObjectGetComm((PetscObject)pc,&comm);
1194:   if (jac->n_local > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only a single block in a subcommunicator is supported");
1195:   jac->n_local = 1; /* currently only a single block is supported for a subcommunicator */
1196:   if (!pc->setupcalled) {
1197:     wasSetup  = PETSC_FALSE;
1198:     PetscNewLog(pc,&mpjac);
1199:     jac->data = (void*)mpjac;

1201:     /* initialize datastructure mpjac */
1202:     if (!jac->psubcomm) {
1203:       /* Create default contiguous subcommunicatiors if user does not provide them */
1204:       PetscSubcommCreate(comm,&jac->psubcomm);
1205:       PetscSubcommSetNumber(jac->psubcomm,jac->n);
1206:       PetscSubcommSetType(jac->psubcomm,PETSC_SUBCOMM_CONTIGUOUS);
1207:       PetscLogObjectMemory((PetscObject)pc,sizeof(PetscSubcomm));
1208:     }
1209:     mpjac->psubcomm = jac->psubcomm;
1210:     subcomm         = PetscSubcommChild(mpjac->psubcomm);

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

1215:     /* create a new PC that processors in each subcomm have copy of */
1216:     PetscMalloc1(1,&jac->ksp);
1217:     KSPCreate(subcomm,&jac->ksp[0]);
1218:     KSPSetErrorIfNotConverged(jac->ksp[0],pc->erroriffailure);
1219:     PetscObjectIncrementTabLevel((PetscObject)jac->ksp[0],(PetscObject)pc,1);
1220:     PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->ksp[0]);
1221:     KSPSetOperators(jac->ksp[0],mpjac->submats,mpjac->submats);
1222:     KSPGetPC(jac->ksp[0],&mpjac->pc);

1224:     PCGetOptionsPrefix(pc,&prefix);
1225:     KSPSetOptionsPrefix(jac->ksp[0],prefix);
1226:     KSPAppendOptionsPrefix(jac->ksp[0],"sub_");
1227:     /*
1228:       PetscMPIInt rank,subsize,subrank;
1229:       MPI_Comm_rank(comm,&rank);
1230:       MPI_Comm_size(subcomm,&subsize);
1231:       MPI_Comm_rank(subcomm,&subrank);

1233:       MatGetLocalSize(mpjac->submats,&m,NULL);
1234:       MatGetSize(mpjac->submats,&n,NULL);
1235:       PetscSynchronizedPrintf(comm,"[%d], sub-size %d,sub-rank %d\n",rank,subsize,subrank);
1236:       PetscSynchronizedFlush(comm,PETSC_STDOUT);
1237:     */

1239:     /* create dummy vectors xsub and ysub */
1240:     MatGetLocalSize(mpjac->submats,&m,&n);
1241:     VecCreateMPIWithArray(subcomm,1,n,PETSC_DECIDE,NULL,&mpjac->xsub);
1242:     VecCreateMPIWithArray(subcomm,1,m,PETSC_DECIDE,NULL,&mpjac->ysub);
1243: #if defined(PETSC_HAVE_CUDA)
1244:     PetscObjectTypeCompareAny((PetscObject)mpjac->submats,&is_gpumatrix,MATAIJCUSPARSE,MATMPIAIJCUSPARSE,"");
1245:     if (is_gpumatrix) {
1246:       VecSetType(mpjac->xsub,VECMPICUDA);
1247:       VecSetType(mpjac->ysub,VECMPICUDA);
1248:     }
1249: #endif
1250: #if defined(PETSC_HAVE_VIENNACL)
1251:     PetscObjectTypeCompareAny((PetscObject)mpjac->submats,&is_gpumatrix,MATAIJVIENNACL,MATMPIAIJVIENNACL,"");
1252:     if (is_gpumatrix) {
1253:       VecSetType(mpjac->xsub,VECMPIVIENNACL);
1254:       VecSetType(mpjac->ysub,VECMPIVIENNACL);
1255:     }
1256: #endif
1257:     PetscLogObjectParent((PetscObject)pc,(PetscObject)mpjac->xsub);
1258:     PetscLogObjectParent((PetscObject)pc,(PetscObject)mpjac->ysub);

1260:     pc->ops->reset   = PCReset_BJacobi_Multiproc;
1261:     pc->ops->destroy = PCDestroy_BJacobi_Multiproc;
1262:     pc->ops->apply   = PCApply_BJacobi_Multiproc;
1263:   } else { /* pc->setupcalled */
1264:     subcomm = PetscSubcommChild(mpjac->psubcomm);
1265:     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
1266:       /* destroy old matrix blocks, then get new matrix blocks */
1267:       if (mpjac->submats) {MatDestroy(&mpjac->submats);}
1268:       MatGetMultiProcBlock(pc->pmat,subcomm,MAT_INITIAL_MATRIX,&mpjac->submats);
1269:     } else {
1270:       MatGetMultiProcBlock(pc->pmat,subcomm,MAT_REUSE_MATRIX,&mpjac->submats);
1271:     }
1272:     KSPSetOperators(jac->ksp[0],mpjac->submats,mpjac->submats);
1273:   }

1275:   if (!wasSetup && pc->setfromoptionscalled) {
1276:     KSPSetFromOptions(jac->ksp[0]);
1277:   }
1278:   return(0);
1279: }