Actual source code: bjacobi.c

petsc-dev 2014-07-25
Report Typos and Errors
  2: /*
  3:    Defines a block Jacobi preconditioner.
  4: */
  5: #include <petsc-private/pcimpl.h>              /*I "petscpc.h" I*/
  6: #include <../src/ksp/pc/impls/bjacobi/bjacobi.h>

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

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

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

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

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

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

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

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

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

153:   PetscFree(jac->g_lens);
154:   PetscFree(jac->l_lens);
155:   PetscFree(pc->data);
156:   return(0);
157: }


162: static PetscErrorCode PCSetFromOptions_BJacobi(PC pc)
163: {
164:   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;
166:   PetscInt       blocks,i;
167:   PetscBool      flg;

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

186: #include <petscdraw.h>
189: static PetscErrorCode PCView_BJacobi(PC pc,PetscViewer viewer)
190: {
191:   PC_BJacobi           *jac   = (PC_BJacobi*)pc->data;
192:   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
193:   PetscErrorCode       ierr;
194:   PetscMPIInt          rank;
195:   PetscInt             i;
196:   PetscBool            iascii,isstring,isdraw;
197:   PetscViewer          sviewer;

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

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

267: /* -------------------------------------------------------------------------------------*/

271: static PetscErrorCode  PCBJacobiGetSubKSP_BJacobi(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp)
272: {
273:   PC_BJacobi *jac = (PC_BJacobi*)pc->data;;

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

278:   if (n_local) *n_local = jac->n_local;
279:   if (first_local) *first_local = jac->first_local;
280:   *ksp                   = jac->ksp;
281:   jac->same_local_solves = PETSC_FALSE;        /* Assume that local solves are now different;
282:                                                   not necessarily true though!  This flag is
283:                                                   used only for PCView_BJacobi() */
284:   return(0);
285: }

289: static PetscErrorCode  PCBJacobiSetTotalBlocks_BJacobi(PC pc,PetscInt blocks,PetscInt *lens)
290: {
291:   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;

295:   if (pc->setupcalled > 0 && jac->n!=blocks) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Cannot alter number of blocks after PCSetUp()/KSPSetUp() has been called");
296:   jac->n = blocks;
297:   if (!lens) jac->g_lens = 0;
298:   else {
299:     PetscMalloc1(blocks,&jac->g_lens);
300:     PetscLogObjectMemory((PetscObject)pc,blocks*sizeof(PetscInt));
301:     PetscMemcpy(jac->g_lens,lens,blocks*sizeof(PetscInt));
302:   }
303:   return(0);
304: }

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

313:   *blocks = jac->n;
314:   if (lens) *lens = jac->g_lens;
315:   return(0);
316: }

320: static PetscErrorCode  PCBJacobiSetLocalBlocks_BJacobi(PC pc,PetscInt blocks,const PetscInt lens[])
321: {
322:   PC_BJacobi     *jac;

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

328:   jac->n_local = blocks;
329:   if (!lens) jac->l_lens = 0;
330:   else {
331:     PetscMalloc1(blocks,&jac->l_lens);
332:     PetscLogObjectMemory((PetscObject)pc,blocks*sizeof(PetscInt));
333:     PetscMemcpy(jac->l_lens,lens,blocks*sizeof(PetscInt));
334:   }
335:   return(0);
336: }

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

345:   *blocks = jac->n_local;
346:   if (lens) *lens = jac->l_lens;
347:   return(0);
348: }

350: /* -------------------------------------------------------------------------------------*/

354: /*@C
355:    PCBJacobiGetSubKSP - Gets the local KSP contexts for all blocks on
356:    this processor.

358:    Note Collective

360:    Input Parameter:
361: .  pc - the preconditioner context

363:    Output Parameters:
364: +  n_local - the number of blocks on this processor, or NULL
365: .  first_local - the global number of the first block on this processor, or NULL
366: -  ksp - the array of KSP contexts

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

371:    Currently for some matrix implementations only 1 block per processor
372:    is supported.

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

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

380:    Level: advanced

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

384: .seealso: PCBJacobiGetSubKSP()
385: @*/
386: PetscErrorCode  PCBJacobiGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
387: {

392:   PetscUseMethod(pc,"PCBJacobiGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));
393:   return(0);
394: }

398: /*@
399:    PCBJacobiSetTotalBlocks - Sets the global number of blocks for the block
400:    Jacobi preconditioner.

402:    Collective on PC

404:    Input Parameters:
405: +  pc - the preconditioner context
406: .  blocks - the number of blocks
407: -  lens - [optional] integer array containing the size of each block

409:    Options Database Key:
410: .  -pc_bjacobi_blocks <blocks> - Sets the number of global blocks

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

416:    Level: intermediate

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

420: .seealso: PCSetUseAmat(), PCBJacobiSetLocalBlocks()
421: @*/
422: PetscErrorCode  PCBJacobiSetTotalBlocks(PC pc,PetscInt blocks,const PetscInt lens[])
423: {

428:   if (blocks <= 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Must have positive blocks");
429:   PetscTryMethod(pc,"PCBJacobiSetTotalBlocks_C",(PC,PetscInt,const PetscInt[]),(pc,blocks,lens));
430:   return(0);
431: }

435: /*@C
436:    PCBJacobiGetTotalBlocks - Gets the global number of blocks for the block
437:    Jacobi preconditioner.

439:    Not Collective

441:    Input Parameter:
442: .  pc - the preconditioner context

444:    Output parameters:
445: +  blocks - the number of blocks
446: -  lens - integer array containing the size of each block

448:    Level: intermediate

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

452: .seealso: PCSetUseAmat(), PCBJacobiGetLocalBlocks()
453: @*/
454: PetscErrorCode  PCBJacobiGetTotalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
455: {

461:   PetscUseMethod(pc,"PCBJacobiGetTotalBlocks_C",(PC,PetscInt*, const PetscInt *[]),(pc,blocks,lens));
462:   return(0);
463: }

467: /*@
468:    PCBJacobiSetLocalBlocks - Sets the local number of blocks for the block
469:    Jacobi preconditioner.

471:    Not Collective

473:    Input Parameters:
474: +  pc - the preconditioner context
475: .  blocks - the number of blocks
476: -  lens - [optional] integer array containing size of each block

478:    Note:
479:    Currently only a limited number of blocking configurations are supported.

481:    Level: intermediate

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

485: .seealso: PCSetUseAmat(), PCBJacobiSetTotalBlocks()
486: @*/
487: PetscErrorCode  PCBJacobiSetLocalBlocks(PC pc,PetscInt blocks,const PetscInt lens[])
488: {

493:   if (blocks < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have nonegative blocks");
494:   PetscTryMethod(pc,"PCBJacobiSetLocalBlocks_C",(PC,PetscInt,const PetscInt []),(pc,blocks,lens));
495:   return(0);
496: }

500: /*@C
501:    PCBJacobiGetLocalBlocks - Gets the local number of blocks for the block
502:    Jacobi preconditioner.

504:    Not Collective

506:    Input Parameters:
507: +  pc - the preconditioner context
508: .  blocks - the number of blocks
509: -  lens - [optional] integer array containing size of each block

511:    Note:
512:    Currently only a limited number of blocking configurations are supported.

514:    Level: intermediate

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

518: .seealso: PCSetUseAmat(), PCBJacobiGetTotalBlocks()
519: @*/
520: PetscErrorCode  PCBJacobiGetLocalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
521: {

527:   PetscUseMethod(pc,"PCBJacobiGetLocalBlocks_C",(PC,PetscInt*, const PetscInt *[]),(pc,blocks,lens));
528:   return(0);
529: }

531: /* -----------------------------------------------------------------------------------*/

533: /*MC
534:    PCBJACOBI - Use block Jacobi preconditioning, each block is (approximately) solved with
535:            its own KSP object.

537:    Options Database Keys:
538: .  -pc_use_amat - use Amat to apply block of operator in inner Krylov method

540:    Notes: Each processor can have one or more blocks, but a block cannot be shared by more
541:      than one processor. Defaults to one block per processor.

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

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

550:    Level: beginner

552:    Concepts: block Jacobi

554:    Developer Notes: This preconditioner does not currently work with CUDA/CUSP for a couple of reasons.
555:        (1) It creates seq vectors as work vectors that should be cusp
556:        (2) The use of VecPlaceArray() is not handled properly by CUSP (that is it will not know where
557:            the ownership of the vector is so may use wrong values) even if it did know the ownership
558:            it may induce extra copy ups and downs. Satish suggests a VecTransplantArray() to handle two
559:            vectors sharing the same pointer and handling the CUSP side as well instead of VecGetArray()/VecPlaceArray().


562: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
563:            PCASM, PCSetUseAmat(), PCGetUseAmat(), PCBJacobiGetSubKSP(), PCBJacobiSetTotalBlocks(),
564:            PCBJacobiSetLocalBlocks(), PCSetModifySubmatrices()
565: M*/

569: PETSC_EXTERN PetscErrorCode PCCreate_BJacobi(PC pc)
570: {
572:   PetscMPIInt    rank;
573:   PC_BJacobi     *jac;

576:   PetscNewLog(pc,&jac);
577:   MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);

579:   pc->ops->apply           = 0;
580:   pc->ops->applytranspose  = 0;
581:   pc->ops->setup           = PCSetUp_BJacobi;
582:   pc->ops->destroy         = PCDestroy_BJacobi;
583:   pc->ops->setfromoptions  = PCSetFromOptions_BJacobi;
584:   pc->ops->view            = PCView_BJacobi;
585:   pc->ops->applyrichardson = 0;

587:   pc->data               = (void*)jac;
588:   jac->n                 = -1;
589:   jac->n_local           = -1;
590:   jac->first_local       = rank;
591:   jac->ksp               = 0;
592:   jac->same_local_solves = PETSC_TRUE;
593:   jac->g_lens            = 0;
594:   jac->l_lens            = 0;
595:   jac->psubcomm          = 0;

597:   PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiGetSubKSP_C",PCBJacobiGetSubKSP_BJacobi);
598:   PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiSetTotalBlocks_C",PCBJacobiSetTotalBlocks_BJacobi);
599:   PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiGetTotalBlocks_C",PCBJacobiGetTotalBlocks_BJacobi);
600:   PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiSetLocalBlocks_C",PCBJacobiSetLocalBlocks_BJacobi);
601:   PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiGetLocalBlocks_C",PCBJacobiGetLocalBlocks_BJacobi);
602:   return(0);
603: }

605: /* --------------------------------------------------------------------------------------------*/
606: /*
607:         These are for a single block per processor; works for AIJ, BAIJ; Seq and MPI
608: */
611: PetscErrorCode PCReset_BJacobi_Singleblock(PC pc)
612: {
613:   PC_BJacobi             *jac  = (PC_BJacobi*)pc->data;
614:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
615:   PetscErrorCode         ierr;

618:   KSPReset(jac->ksp[0]);
619:   VecDestroy(&bjac->x);
620:   VecDestroy(&bjac->y);
621:   return(0);
622: }

626: PetscErrorCode PCDestroy_BJacobi_Singleblock(PC pc)
627: {
628:   PC_BJacobi             *jac  = (PC_BJacobi*)pc->data;
629:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
630:   PetscErrorCode         ierr;

633:   PCReset_BJacobi_Singleblock(pc);
634:   KSPDestroy(&jac->ksp[0]);
635:   PetscFree(jac->ksp);
636:   PetscFree(jac->l_lens);
637:   PetscFree(jac->g_lens);
638:   PetscFree(bjac);
639:   PetscFree(pc->data);
640:   return(0);
641: }

645: PetscErrorCode PCSetUpOnBlocks_BJacobi_Singleblock(PC pc)
646: {
648:   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;

651:   KSPSetUp(jac->ksp[0]);
652:   return(0);
653: }

657: PetscErrorCode PCApply_BJacobi_Singleblock(PC pc,Vec x,Vec y)
658: {
659:   PetscErrorCode         ierr;
660:   PC_BJacobi             *jac  = (PC_BJacobi*)pc->data;
661:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
662:   PetscScalar            *x_array,*y_array;

665:   /*
666:       The VecPlaceArray() is to avoid having to copy the
667:     y vector into the bjac->x vector. The reason for
668:     the bjac->x vector is that we need a sequential vector
669:     for the sequential solve.
670:   */
671:   VecGetArray(x,&x_array);
672:   VecGetArray(y,&y_array);
673:   VecPlaceArray(bjac->x,x_array);
674:   VecPlaceArray(bjac->y,y_array);
675:   KSPSolve(jac->ksp[0],bjac->x,bjac->y);
676:   VecResetArray(bjac->x);
677:   VecResetArray(bjac->y);
678:   VecRestoreArray(x,&x_array);
679:   VecRestoreArray(y,&y_array);
680:   return(0);
681: }

685: PetscErrorCode PCApplySymmetricLeft_BJacobi_Singleblock(PC pc,Vec x,Vec y)
686: {
687:   PetscErrorCode         ierr;
688:   PC_BJacobi             *jac  = (PC_BJacobi*)pc->data;
689:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
690:   PetscScalar            *x_array,*y_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:   VecGetArray(x,&x_array);
701:   VecGetArray(y,&y_array);
702:   VecPlaceArray(bjac->x,x_array);
703:   VecPlaceArray(bjac->y,y_array);
704:   /* apply the symmetric left portion of the inner PC operator */
705:   /* note this by-passes the inner KSP and its options completely */
706:   KSPGetPC(jac->ksp[0],&subpc);
707:   PCApplySymmetricLeft(subpc,bjac->x,bjac->y);
708:   VecResetArray(bjac->x);
709:   VecResetArray(bjac->y);
710:   VecRestoreArray(x,&x_array);
711:   VecRestoreArray(y,&y_array);
712:   return(0);
713: }

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

726:   /*
727:       The VecPlaceArray() is to avoid having to copy the
728:     y vector into the bjac->x vector. The reason for
729:     the bjac->x vector is that we need a sequential vector
730:     for the sequential solve.
731:   */
732:   VecGetArray(x,&x_array);
733:   VecGetArray(y,&y_array);
734:   VecPlaceArray(bjac->x,x_array);
735:   VecPlaceArray(bjac->y,y_array);

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

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

743:   VecRestoreArray(x,&x_array);
744:   VecRestoreArray(y,&y_array);
745:   return(0);
746: }

750: PetscErrorCode PCApplyTranspose_BJacobi_Singleblock(PC pc,Vec x,Vec y)
751: {
752:   PetscErrorCode         ierr;
753:   PC_BJacobi             *jac  = (PC_BJacobi*)pc->data;
754:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
755:   PetscScalar            *x_array,*y_array;

758:   /*
759:       The VecPlaceArray() is to avoid having to copy the
760:     y vector into the bjac->x vector. The reason for
761:     the bjac->x vector is that we need a sequential vector
762:     for the sequential solve.
763:   */
764:   VecGetArray(x,&x_array);
765:   VecGetArray(y,&y_array);
766:   VecPlaceArray(bjac->x,x_array);
767:   VecPlaceArray(bjac->y,y_array);
768:   KSPSolveTranspose(jac->ksp[0],bjac->x,bjac->y);
769:   VecResetArray(bjac->x);
770:   VecResetArray(bjac->y);
771:   VecRestoreArray(x,&x_array);
772:   VecRestoreArray(y,&y_array);
773:   return(0);
774: }

778: static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC pc,Mat mat,Mat pmat)
779: {
780:   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
781:   PetscErrorCode         ierr;
782:   PetscInt               m;
783:   KSP                    ksp;
784:   PC_BJacobi_Singleblock *bjac;
785:   PetscBool              wasSetup = PETSC_TRUE;

788:   if (!pc->setupcalled) {
789:     const char *prefix;

791:     if (!jac->ksp) {
792:       wasSetup = PETSC_FALSE;

794:       KSPCreate(PETSC_COMM_SELF,&ksp);
795:       PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);
796:       PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);
797:       KSPSetType(ksp,KSPPREONLY);
798:       PCGetOptionsPrefix(pc,&prefix);
799:       KSPSetOptionsPrefix(ksp,prefix);
800:       KSPAppendOptionsPrefix(ksp,"sub_");

802:       pc->ops->reset               = PCReset_BJacobi_Singleblock;
803:       pc->ops->destroy             = PCDestroy_BJacobi_Singleblock;
804:       pc->ops->apply               = PCApply_BJacobi_Singleblock;
805:       pc->ops->applysymmetricleft  = PCApplySymmetricLeft_BJacobi_Singleblock;
806:       pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Singleblock;
807:       pc->ops->applytranspose      = PCApplyTranspose_BJacobi_Singleblock;
808:       pc->ops->setuponblocks       = PCSetUpOnBlocks_BJacobi_Singleblock;

810:       PetscMalloc(sizeof(KSP),&jac->ksp);
811:       jac->ksp[0] = ksp;

813:       PetscNewLog(pc,&bjac);
814:       jac->data = (void*)bjac;
815:     } else {
816:       ksp  = jac->ksp[0];
817:       bjac = (PC_BJacobi_Singleblock*)jac->data;
818:     }

820:     /*
821:       The reason we need to generate these vectors is to serve
822:       as the right-hand side and solution vector for the solve on the
823:       block. We do not need to allocate space for the vectors since
824:       that is provided via VecPlaceArray() just before the call to
825:       KSPSolve() on the block.
826:     */
827:     MatGetSize(pmat,&m,&m);
828:     VecCreateSeqWithArray(PETSC_COMM_SELF,1,m,NULL,&bjac->x);
829:     VecCreateSeqWithArray(PETSC_COMM_SELF,1,m,NULL,&bjac->y);
830:     PetscLogObjectParent((PetscObject)pc,(PetscObject)bjac->x);
831:     PetscLogObjectParent((PetscObject)pc,(PetscObject)bjac->y);
832:   } else {
833:     ksp  = jac->ksp[0];
834:     bjac = (PC_BJacobi_Singleblock*)jac->data;
835:   }
836:   if (pc->useAmat) {
837:     KSPSetOperators(ksp,mat,pmat);
838:   } else {
839:     KSPSetOperators(ksp,pmat,pmat);
840:   }
841:   if (!wasSetup && pc->setfromoptionscalled) {
842:     /* If PCSetFromOptions_BJacobi is called later, KSPSetFromOptions will be called at that time. */
843:     KSPSetFromOptions(ksp);
844:   }
845:   return(0);
846: }

848: /* ---------------------------------------------------------------------------------------------*/
851: PetscErrorCode PCReset_BJacobi_Multiblock(PC pc)
852: {
853:   PC_BJacobi            *jac  = (PC_BJacobi*)pc->data;
854:   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
855:   PetscErrorCode        ierr;
856:   PetscInt              i;

859:   if (bjac && bjac->pmat) {
860:     MatDestroyMatrices(jac->n_local,&bjac->pmat);
861:     if (pc->useAmat) {
862:       MatDestroyMatrices(jac->n_local,&bjac->mat);
863:     }
864:   }

866:   for (i=0; i<jac->n_local; i++) {
867:     KSPReset(jac->ksp[i]);
868:     if (bjac && bjac->x) {
869:       VecDestroy(&bjac->x[i]);
870:       VecDestroy(&bjac->y[i]);
871:       ISDestroy(&bjac->is[i]);
872:     }
873:   }
874:   PetscFree(jac->l_lens);
875:   PetscFree(jac->g_lens);
876:   return(0);
877: }

881: PetscErrorCode PCDestroy_BJacobi_Multiblock(PC pc)
882: {
883:   PC_BJacobi            *jac  = (PC_BJacobi*)pc->data;
884:   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
885:   PetscErrorCode        ierr;
886:   PetscInt              i;

889:   PCReset_BJacobi_Multiblock(pc);
890:   if (bjac) {
891:     PetscFree2(bjac->x,bjac->y);
892:     PetscFree(bjac->starts);
893:     PetscFree(bjac->is);
894:   }
895:   PetscFree(jac->data);
896:   for (i=0; i<jac->n_local; i++) {
897:     KSPDestroy(&jac->ksp[i]);
898:   }
899:   PetscFree(jac->ksp);
900:   PetscFree(pc->data);
901:   return(0);
902: }

906: PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiblock(PC pc)
907: {
908:   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;
910:   PetscInt       i,n_local = jac->n_local;

913:   for (i=0; i<n_local; i++) {
914:     KSPSetUp(jac->ksp[i]);
915:   }
916:   return(0);
917: }

919: /*
920:       Preconditioner for block Jacobi
921: */
924: PetscErrorCode PCApply_BJacobi_Multiblock(PC pc,Vec x,Vec y)
925: {
926:   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
927:   PetscErrorCode        ierr;
928:   PetscInt              i,n_local = jac->n_local;
929:   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
930:   PetscScalar           *xin,*yin;

933:   VecGetArray(x,&xin);
934:   VecGetArray(y,&yin);
935:   for (i=0; i<n_local; i++) {
936:     /*
937:        To avoid copying the subvector from x into a workspace we instead
938:        make the workspace vector array point to the subpart of the array of
939:        the global vector.
940:     */
941:     VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);
942:     VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);

944:     PetscLogEventBegin(PC_ApplyOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);
945:     KSPSolve(jac->ksp[i],bjac->x[i],bjac->y[i]);
946:     PetscLogEventEnd(PC_ApplyOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);

948:     VecResetArray(bjac->x[i]);
949:     VecResetArray(bjac->y[i]);
950:   }
951:   VecRestoreArray(x,&xin);
952:   VecRestoreArray(y,&yin);
953:   return(0);
954: }

956: /*
957:       Preconditioner for block Jacobi
958: */
961: PetscErrorCode PCApplyTranspose_BJacobi_Multiblock(PC pc,Vec x,Vec y)
962: {
963:   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
964:   PetscErrorCode        ierr;
965:   PetscInt              i,n_local = jac->n_local;
966:   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
967:   PetscScalar           *xin,*yin;

970:   VecGetArray(x,&xin);
971:   VecGetArray(y,&yin);
972:   for (i=0; i<n_local; i++) {
973:     /*
974:        To avoid copying the subvector from x into a workspace we instead
975:        make the workspace vector array point to the subpart of the array of
976:        the global vector.
977:     */
978:     VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);
979:     VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);

981:     PetscLogEventBegin(PC_ApplyTransposeOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);
982:     KSPSolveTranspose(jac->ksp[i],bjac->x[i],bjac->y[i]);
983:     PetscLogEventEnd(PC_ApplyTransposeOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);

985:     VecResetArray(bjac->x[i]);
986:     VecResetArray(bjac->y[i]);
987:   }
988:   VecRestoreArray(x,&xin);
989:   VecRestoreArray(y,&yin);
990:   return(0);
991: }

995: static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC pc,Mat mat,Mat pmat)
996: {
997:   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
998:   PetscErrorCode        ierr;
999:   PetscInt              m,n_local,N,M,start,i;
1000:   const char            *prefix,*pprefix,*mprefix;
1001:   KSP                   ksp;
1002:   Vec                   x,y;
1003:   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
1004:   PC                    subpc;
1005:   IS                    is;
1006:   MatReuse              scall;

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

1011:   n_local = jac->n_local;

1013:   if (pc->useAmat) {
1014:     PetscBool same;
1015:     PetscObjectTypeCompare((PetscObject)mat,((PetscObject)pmat)->type_name,&same);
1016:     if (!same) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Matrices not of same type");
1017:   }

1019:   if (!pc->setupcalled) {
1020:     scall = MAT_INITIAL_MATRIX;

1022:     if (!jac->ksp) {
1023:       pc->ops->reset         = PCReset_BJacobi_Multiblock;
1024:       pc->ops->destroy       = PCDestroy_BJacobi_Multiblock;
1025:       pc->ops->apply         = PCApply_BJacobi_Multiblock;
1026:       pc->ops->applytranspose= PCApplyTranspose_BJacobi_Multiblock;
1027:       pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiblock;

1029:       PetscNewLog(pc,&bjac);
1030:       PetscMalloc1(n_local,&jac->ksp);
1031:       PetscLogObjectMemory((PetscObject)pc,sizeof(n_local*sizeof(KSP)));
1032:       PetscMalloc2(n_local,&bjac->x,n_local,&bjac->y);
1033:       PetscMalloc1(n_local,&bjac->starts);
1034:       PetscLogObjectMemory((PetscObject)pc,sizeof(n_local*sizeof(PetscScalar)));

1036:       jac->data = (void*)bjac;
1037:       PetscMalloc1(n_local,&bjac->is);
1038:       PetscLogObjectMemory((PetscObject)pc,sizeof(n_local*sizeof(IS)));

1040:       for (i=0; i<n_local; i++) {
1041:         KSPCreate(PETSC_COMM_SELF,&ksp);
1042:         PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);
1043:         PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);
1044:         KSPSetType(ksp,KSPPREONLY);
1045:         KSPGetPC(ksp,&subpc);
1046:         PCGetOptionsPrefix(pc,&prefix);
1047:         KSPSetOptionsPrefix(ksp,prefix);
1048:         KSPAppendOptionsPrefix(ksp,"sub_");

1050:         jac->ksp[i] = ksp;
1051:       }
1052:     } else {
1053:       bjac = (PC_BJacobi_Multiblock*)jac->data;
1054:     }

1056:     start = 0;
1057:     for (i=0; i<n_local; i++) {
1058:       m = jac->l_lens[i];
1059:       /*
1060:       The reason we need to generate these vectors is to serve
1061:       as the right-hand side and solution vector for the solve on the
1062:       block. We do not need to allocate space for the vectors since
1063:       that is provided via VecPlaceArray() just before the call to
1064:       KSPSolve() on the block.

1066:       */
1067:       VecCreateSeq(PETSC_COMM_SELF,m,&x);
1068:       VecCreateSeqWithArray(PETSC_COMM_SELF,1,m,NULL,&y);
1069:       PetscLogObjectParent((PetscObject)pc,(PetscObject)x);
1070:       PetscLogObjectParent((PetscObject)pc,(PetscObject)y);

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

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

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

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

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

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

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

1145: static PetscErrorCode PCDestroy_BJacobi_Multiproc(PC pc)
1146: {
1147:   PC_BJacobi           *jac   = (PC_BJacobi*)pc->data;
1148:   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1149:   PetscErrorCode       ierr;

1152:   PCReset_BJacobi_Multiproc(pc);
1153:   KSPDestroy(&jac->ksp[0]);
1154:   PetscFree(jac->ksp);
1155:   PetscSubcommDestroy(&mpjac->psubcomm);

1157:   PetscFree(mpjac);
1158:   PetscFree(pc->data);
1159:   return(0);
1160: }

1164: static PetscErrorCode PCApply_BJacobi_Multiproc(PC pc,Vec x,Vec y)
1165: {
1166:   PC_BJacobi           *jac   = (PC_BJacobi*)pc->data;
1167:   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1168:   PetscErrorCode       ierr;
1169:   PetscScalar          *xarray,*yarray;

1172:   /* place x's and y's local arrays into xsub and ysub */
1173:   VecGetArray(x,&xarray);
1174:   VecGetArray(y,&yarray);
1175:   VecPlaceArray(mpjac->xsub,xarray);
1176:   VecPlaceArray(mpjac->ysub,yarray);

1178:   /* apply preconditioner on each matrix block */
1179:   PetscLogEventBegin(PC_ApplyOnMproc,jac->ksp[0],mpjac->xsub,mpjac->ysub,0);
1180:   KSPSolve(jac->ksp[0],mpjac->xsub,mpjac->ysub);
1181:   PetscLogEventEnd(PC_ApplyOnMproc,jac->ksp[0],mpjac->xsub,mpjac->ysub,0);

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

1190: #include <petsc-private/matimpl.h>
1193: static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC pc)
1194: {
1195:   PC_BJacobi           *jac   = (PC_BJacobi*)pc->data;
1196:   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1197:   PetscErrorCode       ierr;
1198:   PetscInt             m,n;
1199:   MPI_Comm             comm,subcomm=0;
1200:   const char           *prefix;
1201:   PetscBool            wasSetup = PETSC_TRUE;

1204:   PetscObjectGetComm((PetscObject)pc,&comm);
1205:   if (jac->n_local > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only a single block in a subcommunicator is supported");
1206:   jac->n_local = 1; /* currently only a single block is supported for a subcommunicator */
1207:   if (!pc->setupcalled) {
1208:     wasSetup  = PETSC_FALSE;
1209:     PetscNewLog(pc,&mpjac);
1210:     jac->data = (void*)mpjac;

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

1223:     /* Get matrix blocks of pmat */
1224:     if (!pc->pmat->ops->getmultiprocblock) SETERRQ(PetscObjectComm((PetscObject)pc->pmat),PETSC_ERR_SUP,"No support for the requested operation");
1225:     (*pc->pmat->ops->getmultiprocblock)(pc->pmat,subcomm,MAT_INITIAL_MATRIX,&mpjac->submats);

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

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

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

1250:     /* create dummy vectors xsub and ysub */
1251:     MatGetLocalSize(mpjac->submats,&m,&n);
1252:     VecCreateMPIWithArray(subcomm,1,n,PETSC_DECIDE,NULL,&mpjac->xsub);
1253:     VecCreateMPIWithArray(subcomm,1,m,PETSC_DECIDE,NULL,&mpjac->ysub);
1254:     PetscLogObjectParent((PetscObject)pc,(PetscObject)mpjac->xsub);
1255:     PetscLogObjectParent((PetscObject)pc,(PetscObject)mpjac->ysub);

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

1272:   if (!wasSetup && pc->setfromoptionscalled) {
1273:     KSPSetFromOptions(jac->ksp[0]);
1274:   }
1275:   return(0);
1276: }