Actual source code: bjacobi.c

  2: /*
  3:    Defines a block Jacobi preconditioner.
  4: */
  5: #include <private/pcimpl.h>              /*I "petscpc.h" I*/
  6: #include <../src/ksp/pc/impls/bjacobi/bjacobi.h>

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

114:   /* -------------------------
115:       Determines mat and pmat 
116:   ---------------------------*/
117:   PetscObjectQueryFunction((PetscObject)pc->mat,"MatGetDiagonalBlock_C",(void (**)(void))&f);
118:   if (!f && size == 1) {
119:     mat  = pc->mat;
120:     pmat = pc->pmat;
121:   } else {
122:     if (jac->use_true_local) {
123:       /* use block from true matrix, not preconditioner matrix for local MatMult() */
124:       MatGetDiagonalBlock(pc->mat,&mat);
125:       /* make submatrix have same prefix as entire matrix */
126:       PetscObjectGetOptionsPrefix((PetscObject)pc->mat,&mprefix);
127:       PetscObjectSetOptionsPrefix((PetscObject)mat,mprefix);
128:     }
129:     if (pc->pmat != pc->mat || !jac->use_true_local) {
130:       MatGetDiagonalBlock(pc->pmat,&pmat);
131:       /* make submatrix have same prefix as entire matrix */
132:       PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);
133:       PetscObjectSetOptionsPrefix((PetscObject)pmat,pprefix);
134:     } else {
135:       pmat = mat;
136:     }
137:   }

139:   /* ------
140:      Setup code depends on the number of blocks
141:   */
142:   if (jac->n_local == 0) {
143:     PCSetUp_BJacobi_Singleblock(pc,mat,pmat);
144:   } else {
145:     PCSetUp_BJacobi_Multiblock(pc,mat,pmat);
146:   }
147:   return(0);
148: }

150: /* Default destroy, if it has never been setup */
153: static PetscErrorCode PCDestroy_BJacobi(PC pc)
154: {
155:   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;

159:   PetscFree(jac->g_lens);
160:   PetscFree(jac->l_lens);
161:   PetscFree(pc->data);
162:   return(0);
163: }


168: static PetscErrorCode PCSetFromOptions_BJacobi(PC pc)
169: {
170:   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;
172:   PetscInt       blocks;
173:   PetscBool      flg;

176:   PetscOptionsHead("Block Jacobi options");
177:     PetscOptionsInt("-pc_bjacobi_blocks","Total number of blocks","PCBJacobiSetTotalBlocks",jac->n,&blocks,&flg);
178:     if (flg) {
179:       PCBJacobiSetTotalBlocks(pc,blocks,PETSC_NULL);
180:     }
181:     flg  = PETSC_FALSE;
182:     PetscOptionsBool("-pc_bjacobi_truelocal","Use the true matrix, not preconditioner matrix to define matrix vector product in sub-problems","PCBJacobiSetUseTrueLocal",flg,&flg,PETSC_NULL);
183:     if (flg) {
184:       PCBJacobiSetUseTrueLocal(pc);
185:     }
186:   PetscOptionsTail();
187:   return(0);
188: }

192: static PetscErrorCode PCView_BJacobi(PC pc,PetscViewer viewer)
193: {
194:   PC_BJacobi           *jac = (PC_BJacobi*)pc->data;
195:   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
196:   PetscErrorCode       ierr;
197:   PetscMPIInt          rank;
198:   PetscInt             i;
199:   PetscBool            iascii,isstring;
200:   PetscViewer          sviewer;

203:   PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
204:   PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
205:   if (iascii) {
206:     if (jac->use_true_local) {
207:       PetscViewerASCIIPrintf(viewer,"  block Jacobi: using true local matrix, number of blocks = %D\n",jac->n);
208:     }
209:     PetscViewerASCIIPrintf(viewer,"  block Jacobi: number of blocks = %D\n",jac->n);
210:     MPI_Comm_rank(((PetscObject)pc)->comm,&rank);
211:     if (jac->same_local_solves) {
212:       PetscViewerASCIIPrintf(viewer,"  Local solve is same for all blocks, in the following KSP and PC objects:\n");
213:       if (jac->ksp) {
214:         PetscViewerGetSingleton(viewer,&sviewer);
215:         if (!rank){
216:           PetscViewerASCIIPushTab(viewer);
217:           KSPView(jac->ksp[0],sviewer);
218:           PetscViewerASCIIPopTab(viewer);
219:         }
220:         PetscViewerRestoreSingleton(viewer,&sviewer);
221:       } else if (jac->psubcomm && !jac->psubcomm->color){
222:         PetscViewerASCIIGetStdout(mpjac->psubcomm->comm,&sviewer);
223:         PetscViewerASCIIPushTab(viewer);
224:         KSPView(mpjac->ksp,sviewer);
225:         PetscViewerASCIIPopTab(viewer);
226:       }
227:     } else {
228:       PetscInt n_global;
229:       MPI_Allreduce(&jac->n_local,&n_global,1,MPIU_INT,MPI_MAX,((PetscObject)pc)->comm);
230:       PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);
231:       PetscViewerASCIIPrintf(viewer,"  Local solve info for each block is in the following KSP and PC objects:\n");
232:       PetscViewerASCIISynchronizedPrintf(viewer,"[%d] number of local blocks = %D, first local block number = %D\n",
233:                    rank,jac->n_local,jac->first_local);
234:       PetscViewerASCIIPushTab(viewer);
235:       for (i=0; i<n_global; i++) {
236:         PetscViewerGetSingleton(viewer,&sviewer);
237:         if (i < jac->n_local) {
238:           PetscViewerASCIISynchronizedPrintf(viewer,"[%d] local block number %D\n",rank,i);
239:           KSPView(jac->ksp[i],sviewer);
240:           PetscViewerASCIISynchronizedPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");
241:         }
242:         PetscViewerRestoreSingleton(viewer,&sviewer);
243:       }
244:       PetscViewerASCIIPopTab(viewer);
245:       PetscViewerFlush(viewer);
246:       PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);
247:     }
248:   } else if (isstring) {
249:     PetscViewerStringSPrintf(viewer," blks=%D",jac->n);
250:     PetscViewerGetSingleton(viewer,&sviewer);
251:     if (jac->ksp) {KSPView(jac->ksp[0],sviewer);}
252:     PetscViewerRestoreSingleton(viewer,&sviewer);
253:   } else {
254:     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for block Jacobi",((PetscObject)viewer)->type_name);
255:   }
256:   return(0);
257: }

259: /* -------------------------------------------------------------------------------------*/

264: PetscErrorCode  PCBJacobiSetUseTrueLocal_BJacobi(PC pc)
265: {
266:   PC_BJacobi   *jac;

269:   jac                 = (PC_BJacobi*)pc->data;
270:   jac->use_true_local = PETSC_TRUE;
271:   return(0);
272: }

278: PetscErrorCode  PCBJacobiGetSubKSP_BJacobi(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp)
279: {
280:   PC_BJacobi   *jac = (PC_BJacobi*)pc->data;;

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

285:   if (n_local)     *n_local     = jac->n_local;
286:   if (first_local) *first_local = jac->first_local;
287:   *ksp                          = jac->ksp;
288:   jac->same_local_solves        = PETSC_FALSE; /* Assume that local solves are now different;
289:                                                   not necessarily true though!  This flag is 
290:                                                   used only for PCView_BJacobi() */
291:   return(0);
292: }

298: PetscErrorCode  PCBJacobiSetTotalBlocks_BJacobi(PC pc,PetscInt blocks,PetscInt *lens)
299: {
300:   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;


305:   if (pc->setupcalled > 0 && jac->n!=blocks) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ORDER,"Cannot alter number of blocks after PCSetUp()/KSPSetUp() has been called");
306:   jac->n = blocks;
307:   if (!lens) {
308:     jac->g_lens = 0;
309:   } else {
310:     PetscMalloc(blocks*sizeof(PetscInt),&jac->g_lens);
311:     PetscLogObjectMemory(pc,blocks*sizeof(PetscInt));
312:     PetscMemcpy(jac->g_lens,lens,blocks*sizeof(PetscInt));
313:   }
314:   return(0);
315: }

321: PetscErrorCode  PCBJacobiGetTotalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[])
322: {
323:   PC_BJacobi *jac = (PC_BJacobi*) pc->data;

326:   *blocks = jac->n;
327:   if (lens) *lens = jac->g_lens;
328:   return(0);
329: }

335: PetscErrorCode  PCBJacobiSetLocalBlocks_BJacobi(PC pc,PetscInt blocks,const PetscInt lens[])
336: {
337:   PC_BJacobi     *jac;

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

343:   jac->n_local = blocks;
344:   if (!lens) {
345:     jac->l_lens = 0;
346:   } else {
347:     PetscMalloc(blocks*sizeof(PetscInt),&jac->l_lens);
348:     PetscLogObjectMemory(pc,blocks*sizeof(PetscInt));
349:     PetscMemcpy(jac->l_lens,lens,blocks*sizeof(PetscInt));
350:   }
351:   return(0);
352: }

358: PetscErrorCode  PCBJacobiGetLocalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[])
359: {
360:   PC_BJacobi *jac = (PC_BJacobi*) pc->data;

363:   *blocks = jac->n_local;
364:   if (lens) *lens = jac->l_lens;
365:   return(0);
366: }

369: /* -------------------------------------------------------------------------------------*/

373: /*@
374:    PCBJacobiSetUseTrueLocal - Sets a flag to indicate that the block 
375:    problem is associated with the linear system matrix instead of the
376:    default (where it is associated with the preconditioning matrix).
377:    That is, if the local system is solved iteratively then it iterates
378:    on the block from the matrix using the block from the preconditioner
379:    as the preconditioner for the local block.

381:    Logically Collective on PC

383:    Input Parameters:
384: .  pc - the preconditioner context

386:    Options Database Key:
387: .  -pc_bjacobi_truelocal - Activates PCBJacobiSetUseTrueLocal()

389:    Notes:
390:    For the common case in which the preconditioning and linear 
391:    system matrices are identical, this routine is unnecessary.

393:    Level: intermediate

395: .keywords:  block, Jacobi, set, true, local, flag

397: .seealso: PCSetOperators(), PCBJacobiSetLocalBlocks()
398: @*/
399: PetscErrorCode  PCBJacobiSetUseTrueLocal(PC pc)
400: {

405:   PetscTryMethod(pc,"PCBJacobiSetUseTrueLocal_C",(PC),(pc));
406:   return(0);
407: }

411: /*@C
412:    PCBJacobiGetSubKSP - Gets the local KSP contexts for all blocks on
413:    this processor.
414:    
415:    Note Collective

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

420:    Output Parameters:
421: +  n_local - the number of blocks on this processor, or PETSC_NULL
422: .  first_local - the global number of the first block on this processor, or PETSC_NULL
423: -  ksp - the array of KSP contexts

425:    Notes:  
426:    After PCBJacobiGetSubKSP() the array of KSP contexts is not to be freed.
427:    
428:    Currently for some matrix implementations only 1 block per processor 
429:    is supported.
430:    
431:    You must call KSPSetUp() or PCSetUp() before calling PCBJacobiGetSubKSP().

433:    Level: advanced

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

437: .seealso: PCBJacobiGetSubKSP()
438: @*/
439: PetscErrorCode  PCBJacobiGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
440: {

445:   PetscUseMethod(pc,"PCBJacobiGetSubKSP_C",(PC,PetscInt *,PetscInt *,KSP **),(pc,n_local,first_local,ksp));
446:   return(0);
447: }

451: /*@
452:    PCBJacobiSetTotalBlocks - Sets the global number of blocks for the block
453:    Jacobi preconditioner.

455:    Collective on PC

457:    Input Parameters:
458: +  pc - the preconditioner context
459: .  blocks - the number of blocks
460: -  lens - [optional] integer array containing the size of each block

462:    Options Database Key:
463: .  -pc_bjacobi_blocks <blocks> - Sets the number of global blocks

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

469:    Level: intermediate

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

473: .seealso: PCBJacobiSetUseTrueLocal(), PCBJacobiSetLocalBlocks()
474: @*/
475: PetscErrorCode  PCBJacobiSetTotalBlocks(PC pc,PetscInt blocks,const PetscInt lens[])
476: {

481:   if (blocks <= 0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have positive blocks");
482:   PetscTryMethod(pc,"PCBJacobiSetTotalBlocks_C",(PC,PetscInt,const PetscInt[]),(pc,blocks,lens));
483:   return(0);
484: }

488: /*@C
489:    PCBJacobiGetTotalBlocks - Gets the global number of blocks for the block
490:    Jacobi preconditioner.

492:    Not Collective

494:    Input Parameter:
495: .  pc - the preconditioner context

497:    Output parameters:
498: +  blocks - the number of blocks
499: -  lens - integer array containing the size of each block

501:    Level: intermediate

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

505: .seealso: PCBJacobiSetUseTrueLocal(), PCBJacobiGetLocalBlocks()
506: @*/
507: PetscErrorCode  PCBJacobiGetTotalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
508: {

514:   PetscUseMethod(pc,"PCBJacobiGetTotalBlocks_C",(PC,PetscInt*, const PetscInt *[]),(pc,blocks,lens));
515:   return(0);
516: }
517: 
520: /*@
521:    PCBJacobiSetLocalBlocks - Sets the local number of blocks for the block
522:    Jacobi preconditioner.

524:    Not Collective

526:    Input Parameters:
527: +  pc - the preconditioner context
528: .  blocks - the number of blocks
529: -  lens - [optional] integer array containing size of each block

531:    Note:  
532:    Currently only a limited number of blocking configurations are supported.

534:    Level: intermediate

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

538: .seealso: PCBJacobiSetUseTrueLocal(), PCBJacobiSetTotalBlocks()
539: @*/
540: PetscErrorCode  PCBJacobiSetLocalBlocks(PC pc,PetscInt blocks,const PetscInt lens[])
541: {

546:   if (blocks < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have nonegative blocks");
547:   PetscTryMethod(pc,"PCBJacobiSetLocalBlocks_C",(PC,PetscInt,const PetscInt []),(pc,blocks,lens));
548:   return(0);
549: }
550: 
553: /*@C
554:    PCBJacobiGetLocalBlocks - Gets the local number of blocks for the block
555:    Jacobi preconditioner.

557:    Not Collective

559:    Input Parameters:
560: +  pc - the preconditioner context
561: .  blocks - the number of blocks
562: -  lens - [optional] integer array containing size of each block

564:    Note:  
565:    Currently only a limited number of blocking configurations are supported.

567:    Level: intermediate

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

571: .seealso: PCBJacobiSetUseTrueLocal(), PCBJacobiGetTotalBlocks()
572: @*/
573: PetscErrorCode  PCBJacobiGetLocalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
574: {

580:   PetscUseMethod(pc,"PCBJacobiGetLocalBlocks_C",(PC,PetscInt*, const PetscInt *[]),(pc,blocks,lens));
581:   return(0);
582: }

584: /* -----------------------------------------------------------------------------------*/

586: /*MC
587:    PCBJACOBI - Use block Jacobi preconditioning, each block is (approximately) solved with 
588:            its own KSP object.

590:    Options Database Keys:
591: .  -pc_bjacobi_truelocal - Activates PCBJacobiSetUseTrueLocal()

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

596:      To set options on the solvers for each block append -sub_ to all the KSP, KSP, and PC
597:         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
598:         
599:      To set the options on the solvers separate for each block call PCBJacobiGetSubKSP()
600:          and set the options directly on the resulting KSP object (you can access its PC
601:          KSPGetPC())

603:    Level: beginner

605:    Concepts: block Jacobi

607: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
608:            PCASM, PCBJacobiSetUseTrueLocal(), PCBJacobiGetSubKSP(), PCBJacobiSetTotalBlocks(),
609:            PCBJacobiSetLocalBlocks(), PCSetModifySubmatrices()
610: M*/

615: PetscErrorCode  PCCreate_BJacobi(PC pc)
616: {
618:   PetscMPIInt    rank;
619:   PC_BJacobi     *jac;

622:   PetscNewLog(pc,PC_BJacobi,&jac);
623:   MPI_Comm_rank(((PetscObject)pc)->comm,&rank);
624:   pc->ops->apply              = 0;
625:   pc->ops->applytranspose     = 0;
626:   pc->ops->setup              = PCSetUp_BJacobi;
627:   pc->ops->destroy            = PCDestroy_BJacobi;
628:   pc->ops->setfromoptions     = PCSetFromOptions_BJacobi;
629:   pc->ops->view               = PCView_BJacobi;
630:   pc->ops->applyrichardson    = 0;

632:   pc->data               = (void*)jac;
633:   jac->n                 = -1;
634:   jac->n_local           = -1;
635:   jac->first_local       = rank;
636:   jac->ksp               = 0;
637:   jac->use_true_local    = PETSC_FALSE;
638:   jac->same_local_solves = PETSC_TRUE;
639:   jac->g_lens            = 0;
640:   jac->l_lens            = 0;
641:   jac->psubcomm          = 0;

643:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiSetUseTrueLocal_C",
644:                     "PCBJacobiSetUseTrueLocal_BJacobi",
645:                     PCBJacobiSetUseTrueLocal_BJacobi);
646:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiGetSubKSP_C","PCBJacobiGetSubKSP_BJacobi",
647:                     PCBJacobiGetSubKSP_BJacobi);
648:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiSetTotalBlocks_C","PCBJacobiSetTotalBlocks_BJacobi",
649:                     PCBJacobiSetTotalBlocks_BJacobi);
650:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiGetTotalBlocks_C","PCBJacobiGetTotalBlocks_BJacobi",
651:                     PCBJacobiGetTotalBlocks_BJacobi);
652:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiSetLocalBlocks_C","PCBJacobiSetLocalBlocks_BJacobi",
653:                     PCBJacobiSetLocalBlocks_BJacobi);
654:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiGetLocalBlocks_C","PCBJacobiGetLocalBlocks_BJacobi",
655:                     PCBJacobiGetLocalBlocks_BJacobi);

657:   return(0);
658: }

661: /* --------------------------------------------------------------------------------------------*/
662: /*
663:         These are for a single block per processor; works for AIJ, BAIJ; Seq and MPI
664: */
667: PetscErrorCode PCReset_BJacobi_Singleblock(PC pc)
668: {
669:   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
670:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
671:   PetscErrorCode         ierr;

674:   KSPReset(jac->ksp[0]);
675:   VecDestroy(&bjac->x);
676:   VecDestroy(&bjac->y);
677:   return(0);
678: }

682: PetscErrorCode PCDestroy_BJacobi_Singleblock(PC pc)
683: {
684:   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
685:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
686:   PetscErrorCode         ierr;

689:   PCReset_BJacobi_Singleblock(pc);
690:   KSPDestroy(&jac->ksp[0]);
691:   PetscFree(jac->ksp);
692:   PetscFree(jac->l_lens);
693:   PetscFree(jac->g_lens);
694:   PetscFree(bjac);
695:   PetscFree(pc->data);
696:   return(0);
697: }

701: PetscErrorCode PCSetUpOnBlocks_BJacobi_Singleblock(PC pc)
702: {
704:   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;

707:   KSPSetUp(jac->ksp[0]);
708:   return(0);
709: }

713: PetscErrorCode PCApply_BJacobi_Singleblock(PC pc,Vec x,Vec y)
714: {
715:   PetscErrorCode         ierr;
716:   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
717:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
718:   PetscScalar            *x_array,*y_array;

721:   /* 
722:       The VecPlaceArray() is to avoid having to copy the 
723:     y vector into the bjac->x vector. The reason for 
724:     the bjac->x vector is that we need a sequential vector
725:     for the sequential solve.
726:   */
727:   VecGetArray(x,&x_array);
728:   VecGetArray(y,&y_array);
729:   VecPlaceArray(bjac->x,x_array);
730:   VecPlaceArray(bjac->y,y_array);
731:   KSPSolve(jac->ksp[0],bjac->x,bjac->y);
732:   VecResetArray(bjac->x);
733:   VecResetArray(bjac->y);
734:   VecRestoreArray(x,&x_array);
735:   VecRestoreArray(y,&y_array);
736:   return(0);
737: }

741: PetscErrorCode PCApplySymmetricLeft_BJacobi_Singleblock(PC pc,Vec x,Vec y)
742: {
743:   PetscErrorCode         ierr;
744:   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
745:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
746:   PetscScalar            *x_array,*y_array;
747:   PC                     subpc;

750:   /* 
751:       The VecPlaceArray() is to avoid having to copy the 
752:     y vector into the bjac->x vector. The reason for 
753:     the bjac->x vector is that we need a sequential vector
754:     for the sequential solve.
755:   */
756:   VecGetArray(x,&x_array);
757:   VecGetArray(y,&y_array);
758:   VecPlaceArray(bjac->x,x_array);
759:   VecPlaceArray(bjac->y,y_array);

761:   /* apply the symmetric left portion of the inner PC operator */
762:   /* note this by-passes the inner KSP and its options completely */

764:   KSPGetPC(jac->ksp[0],&subpc);
765:   PCApplySymmetricLeft(subpc,bjac->x,bjac->y);
766:   VecResetArray(bjac->x);
767:   VecResetArray(bjac->y);

769:   VecRestoreArray(x,&x_array);
770:   VecRestoreArray(y,&y_array);
771:   return(0);
772: }

776: PetscErrorCode PCApplySymmetricRight_BJacobi_Singleblock(PC pc,Vec x,Vec y)
777: {
778:   PetscErrorCode         ierr;
779:   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
780:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
781:   PetscScalar            *x_array,*y_array;
782:   PC                     subpc;

785:   /* 
786:       The VecPlaceArray() is to avoid having to copy the 
787:     y vector into the bjac->x vector. The reason for 
788:     the bjac->x vector is that we need a sequential vector
789:     for the sequential solve.
790:   */
791:   VecGetArray(x,&x_array);
792:   VecGetArray(y,&y_array);
793:   VecPlaceArray(bjac->x,x_array);
794:   VecPlaceArray(bjac->y,y_array);

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

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

802:   VecRestoreArray(x,&x_array);
803:   VecRestoreArray(y,&y_array);
804:   return(0);
805: }

809: PetscErrorCode PCApplyTranspose_BJacobi_Singleblock(PC pc,Vec x,Vec y)
810: {
811:   PetscErrorCode         ierr;
812:   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
813:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
814:   PetscScalar            *x_array,*y_array;

817:   /* 
818:       The VecPlaceArray() is to avoid having to copy the 
819:     y vector into the bjac->x vector. The reason for 
820:     the bjac->x vector is that we need a sequential vector
821:     for the sequential solve.
822:   */
823:   VecGetArray(x,&x_array);
824:   VecGetArray(y,&y_array);
825:   VecPlaceArray(bjac->x,x_array);
826:   VecPlaceArray(bjac->y,y_array);
827:   KSPSolveTranspose(jac->ksp[0],bjac->x,bjac->y);
828:   VecResetArray(bjac->x);
829:   VecResetArray(bjac->y);
830:   VecRestoreArray(x,&x_array);
831:   VecRestoreArray(y,&y_array);
832:   return(0);
833: }

837: static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC pc,Mat mat,Mat pmat)
838: {
839:   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
840:   PetscErrorCode         ierr;
841:   PetscInt               m;
842:   KSP                    ksp;
843:   PC_BJacobi_Singleblock *bjac;
844:   PetscBool              wasSetup;


848:   if (!pc->setupcalled) {
849:     const char *prefix;

851:     if (!jac->ksp) {
852:       wasSetup = PETSC_FALSE;
853:       KSPCreate(PETSC_COMM_SELF,&ksp);
854:       PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);
855:       PetscLogObjectParent(pc,ksp);
856:       KSPSetType(ksp,KSPPREONLY);
857:       PCGetOptionsPrefix(pc,&prefix);
858:       KSPSetOptionsPrefix(ksp,prefix);
859:       KSPAppendOptionsPrefix(ksp,"sub_");

861:       pc->ops->reset               = PCReset_BJacobi_Singleblock;
862:       pc->ops->destroy             = PCDestroy_BJacobi_Singleblock;
863:       pc->ops->apply               = PCApply_BJacobi_Singleblock;
864:       pc->ops->applysymmetricleft  = PCApplySymmetricLeft_BJacobi_Singleblock;
865:       pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Singleblock;
866:       pc->ops->applytranspose      = PCApplyTranspose_BJacobi_Singleblock;
867:       pc->ops->setuponblocks       = PCSetUpOnBlocks_BJacobi_Singleblock;

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

872:       PetscNewLog(pc,PC_BJacobi_Singleblock,&bjac);
873:       jac->data    = (void*)bjac;
874:     } else {
875:       ksp  = jac->ksp[0];
876:       bjac = (PC_BJacobi_Singleblock *)jac->data;
877:     }

879:     /*
880:       The reason we need to generate these vectors is to serve 
881:       as the right-hand side and solution vector for the solve on the 
882:       block. We do not need to allocate space for the vectors since
883:       that is provided via VecPlaceArray() just before the call to 
884:       KSPSolve() on the block.
885:     */
886:     MatGetSize(pmat,&m,&m);
887:     VecCreateSeqWithArray(PETSC_COMM_SELF,m,PETSC_NULL,&bjac->x);
888:     VecCreateSeqWithArray(PETSC_COMM_SELF,m,PETSC_NULL,&bjac->y);
889:     PetscLogObjectParent(pc,bjac->x);
890:     PetscLogObjectParent(pc,bjac->y);
891:   } else {
892:     wasSetup = PETSC_TRUE;
893:     ksp = jac->ksp[0];
894:     bjac = (PC_BJacobi_Singleblock *)jac->data;
895:   }
896:   if (jac->use_true_local) {
897:     KSPSetOperators(ksp,mat,pmat,pc->flag);
898:   }  else {
899:     KSPSetOperators(ksp,pmat,pmat,pc->flag);
900:   }
901:   if (!wasSetup && pc->setfromoptionscalled) {
902:     KSPSetFromOptions(ksp);
903:   }
904:   return(0);
905: }

907: /* ---------------------------------------------------------------------------------------------*/
910: PetscErrorCode PCReset_BJacobi_Multiblock(PC pc)
911: {
912:   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
913:   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
914:   PetscErrorCode        ierr;
915:   PetscInt              i;

918:   if (bjac && bjac->pmat) {
919:     MatDestroyMatrices(jac->n_local,&bjac->pmat);
920:     if (jac->use_true_local) {
921:       MatDestroyMatrices(jac->n_local,&bjac->mat);
922:     }
923:   }

925:   for (i=0; i<jac->n_local; i++) {
926:     KSPReset(jac->ksp[i]);
927:     if (bjac && bjac->x) {
928:       VecDestroy(&bjac->x[i]);
929:       VecDestroy(&bjac->y[i]);
930:       ISDestroy(&bjac->is[i]);
931:     }
932:   }
933:   PetscFree(jac->l_lens);
934:   PetscFree(jac->g_lens);
935:   return(0);
936: }

940: PetscErrorCode PCDestroy_BJacobi_Multiblock(PC pc)
941: {
942:   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
943:   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
944:   PetscErrorCode        ierr;
945:   PetscInt              i;

948:   PCReset_BJacobi_Multiblock(pc);
949:   if (bjac) {
950:     PetscFree2(bjac->x,bjac->y);
951:     PetscFree(bjac->starts);
952:     PetscFree(bjac->is);
953:   }
954:   PetscFree(jac->data);
955:   for (i=0; i<jac->n_local; i++) {
956:     KSPDestroy(&jac->ksp[i]);
957:   }
958:   PetscFree(jac->ksp);
959:   PetscFree(pc->data);
960:   return(0);
961: }

965: PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiblock(PC pc)
966: {
967:   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;
969:   PetscInt       i,n_local = jac->n_local;

972:   for (i=0; i<n_local; i++) {
973:     KSPSetUp(jac->ksp[i]);
974:   }
975:   return(0);
976: }

978: /*
979:       Preconditioner for block Jacobi 
980: */
983: PetscErrorCode PCApply_BJacobi_Multiblock(PC pc,Vec x,Vec y)
984: {
985:   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
986:   PetscErrorCode        ierr;
987:   PetscInt              i,n_local = jac->n_local;
988:   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
989:   PetscScalar           *xin,*yin;

992:   VecGetArray(x,&xin);
993:   VecGetArray(y,&yin);
994:   for (i=0; i<n_local; i++) {
995:     /* 
996:        To avoid copying the subvector from x into a workspace we instead 
997:        make the workspace vector array point to the subpart of the array of
998:        the global vector.
999:     */
1000:     VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);
1001:     VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);

1003:     PetscLogEventBegin(PC_SetUpOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);
1004:     KSPSolve(jac->ksp[i],bjac->x[i],bjac->y[i]);
1005:     PetscLogEventEnd(PC_SetUpOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);

1007:     VecResetArray(bjac->x[i]);
1008:     VecResetArray(bjac->y[i]);
1009:   }
1010:   VecRestoreArray(x,&xin);
1011:   VecRestoreArray(y,&yin);
1012:   return(0);
1013: }

1015: /*
1016:       Preconditioner for block Jacobi 
1017: */
1020: PetscErrorCode PCApplyTranspose_BJacobi_Multiblock(PC pc,Vec x,Vec y)
1021: {
1022:   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
1023:   PetscErrorCode        ierr;
1024:   PetscInt              i,n_local = jac->n_local;
1025:   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
1026:   PetscScalar           *xin,*yin;

1029:   VecGetArray(x,&xin);
1030:   VecGetArray(y,&yin);
1031:   for (i=0; i<n_local; i++) {
1032:     /* 
1033:        To avoid copying the subvector from x into a workspace we instead 
1034:        make the workspace vector array point to the subpart of the array of
1035:        the global vector.
1036:     */
1037:     VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);
1038:     VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);

1040:     PetscLogEventBegin(PC_ApplyTransposeOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);
1041:     KSPSolveTranspose(jac->ksp[i],bjac->x[i],bjac->y[i]);
1042:     PetscLogEventEnd(PC_ApplyTransposeOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);

1044:     VecResetArray(bjac->x[i]);
1045:     VecResetArray(bjac->y[i]);
1046:   }
1047:   VecRestoreArray(x,&xin);
1048:   VecRestoreArray(y,&yin);
1049:   return(0);
1050: }

1054: static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC pc,Mat mat,Mat pmat)
1055: {
1056:   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
1057:   PetscErrorCode         ierr;
1058:   PetscInt               m,n_local,N,M,start,i;
1059:   const char             *prefix,*pprefix,*mprefix;
1060:   KSP                    ksp;
1061:   Vec                    x,y;
1062:   PC_BJacobi_Multiblock  *bjac = (PC_BJacobi_Multiblock*)jac->data;
1063:   PC                     subpc;
1064:   IS                     is;
1065:   MatReuse               scall;

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

1070:   n_local = jac->n_local;

1072:   if (jac->use_true_local) {
1073:     PetscBool  same;
1074:     PetscTypeCompare((PetscObject)mat,((PetscObject)pmat)->type_name,&same);
1075:     if (!same) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Matrices not of same type");
1076:   }

1078:   if (!pc->setupcalled) {
1079:     scall                  = MAT_INITIAL_MATRIX;

1081:     if (!jac->ksp) {
1082:       pc->ops->reset         = PCReset_BJacobi_Multiblock;
1083:       pc->ops->destroy       = PCDestroy_BJacobi_Multiblock;
1084:       pc->ops->apply         = PCApply_BJacobi_Multiblock;
1085:       pc->ops->applytranspose= PCApplyTranspose_BJacobi_Multiblock;
1086:       pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiblock;

1088:       PetscNewLog(pc,PC_BJacobi_Multiblock,&bjac);
1089:       PetscMalloc(n_local*sizeof(KSP),&jac->ksp);
1090:       PetscLogObjectMemory(pc,sizeof(n_local*sizeof(KSP)));
1091:       PetscMalloc2(n_local,Vec,&bjac->x,n_local,Vec,&bjac->y);
1092:       PetscMalloc(n_local*sizeof(PetscScalar),&bjac->starts);
1093:       PetscLogObjectMemory(pc,sizeof(n_local*sizeof(PetscScalar)));
1094: 
1095:       jac->data    = (void*)bjac;
1096:       PetscMalloc(n_local*sizeof(IS),&bjac->is);
1097:       PetscLogObjectMemory(pc,sizeof(n_local*sizeof(IS)));

1099:       for (i=0; i<n_local; i++) {
1100:         KSPCreate(PETSC_COMM_SELF,&ksp);
1101:         PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);
1102:         PetscLogObjectParent(pc,ksp);
1103:         KSPSetType(ksp,KSPPREONLY);
1104:         KSPGetPC(ksp,&subpc);
1105:         PCGetOptionsPrefix(pc,&prefix);
1106:         KSPSetOptionsPrefix(ksp,prefix);
1107:         KSPAppendOptionsPrefix(ksp,"sub_");
1108:         jac->ksp[i]    = ksp;
1109:       }
1110:     } else {
1111:       bjac = (PC_BJacobi_Multiblock*)jac->data;
1112:     }

1114:     start = 0;
1115:     for (i=0; i<n_local; i++) {
1116:       m = jac->l_lens[i];
1117:       /*
1118:       The reason we need to generate these vectors is to serve 
1119:       as the right-hand side and solution vector for the solve on the 
1120:       block. We do not need to allocate space for the vectors since
1121:       that is provided via VecPlaceArray() just before the call to 
1122:       KSPSolve() on the block.

1124:       */
1125:       VecCreateSeq(PETSC_COMM_SELF,m,&x);
1126:       VecCreateSeqWithArray(PETSC_COMM_SELF,m,PETSC_NULL,&y);
1127:       PetscLogObjectParent(pc,x);
1128:       PetscLogObjectParent(pc,y);
1129:       bjac->x[i]      = x;
1130:       bjac->y[i]      = y;
1131:       bjac->starts[i] = start;

1133:       ISCreateStride(PETSC_COMM_SELF,m,start,1,&is);
1134:       bjac->is[i] = is;
1135:       PetscLogObjectParent(pc,is);

1137:       start += m;
1138:     }
1139:   } else {
1140:     bjac = (PC_BJacobi_Multiblock*)jac->data;
1141:     /* 
1142:        Destroy the blocks from the previous iteration
1143:     */
1144:     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
1145:       MatDestroyMatrices(n_local,&bjac->pmat);
1146:       if (jac->use_true_local) {
1147:         MatDestroyMatrices(n_local,&bjac->mat);
1148:       }
1149:       scall = MAT_INITIAL_MATRIX;
1150:     } else {
1151:       scall = MAT_REUSE_MATRIX;
1152:     }
1153:   }

1155:   MatGetSubMatrices(pmat,n_local,bjac->is,bjac->is,scall,&bjac->pmat);
1156:   if (jac->use_true_local) {
1157:     PetscObjectGetOptionsPrefix((PetscObject)mat,&mprefix);
1158:     MatGetSubMatrices(mat,n_local,bjac->is,bjac->is,scall,&bjac->mat);
1159:   }
1160:   /* Return control to the user so that the submatrices can be modified (e.g., to apply
1161:      different boundary conditions for the submatrices than for the global problem) */
1162:   PCModifySubMatrices(pc,n_local,bjac->is,bjac->is,bjac->pmat,pc->modifysubmatricesP);

1164:   PetscObjectGetOptionsPrefix((PetscObject)pmat,&pprefix);
1165:   for (i=0; i<n_local; i++) {
1166:     PetscLogObjectParent(pc,bjac->pmat[i]);
1167:     PetscObjectSetOptionsPrefix((PetscObject)bjac->pmat[i],pprefix);
1168:     if (jac->use_true_local) {
1169:       PetscLogObjectParent(pc,bjac->mat[i]);
1170:       PetscObjectSetOptionsPrefix((PetscObject)bjac->mat[i],mprefix);
1171:       KSPSetOperators(jac->ksp[i],bjac->mat[i],bjac->pmat[i],pc->flag);
1172:     } else {
1173:       KSPSetOperators(jac->ksp[i],bjac->pmat[i],bjac->pmat[i],pc->flag);
1174:     }
1175:     if(pc->setfromoptionscalled){
1176:       KSPSetFromOptions(jac->ksp[i]);
1177:     }
1178:   }
1179:   return(0);
1180: }

1182: /* ---------------------------------------------------------------------------------------------*/
1183: /*
1184:       These are for a single block with multiple processes; 
1185: */
1188: static PetscErrorCode PCReset_BJacobi_Multiproc(PC pc)
1189: {
1190:   PC_BJacobi           *jac = (PC_BJacobi*)pc->data;
1191:   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1192:   PetscErrorCode       ierr;

1195:   VecDestroy(&mpjac->ysub);
1196:   VecDestroy(&mpjac->xsub);
1197:   MatDestroy(&mpjac->submats);
1198:   if (mpjac->ksp){KSPReset(mpjac->ksp);}
1199:   return(0);
1200: }

1204: static PetscErrorCode PCDestroy_BJacobi_Multiproc(PC pc)
1205: {
1206:   PC_BJacobi           *jac = (PC_BJacobi*)pc->data;
1207:   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1208:   PetscErrorCode       ierr;

1211:   PCReset_BJacobi_Multiproc(pc);
1212:   KSPDestroy(&mpjac->ksp);
1213:   PetscSubcommDestroy(&mpjac->psubcomm);

1215:   PetscFree(mpjac);
1216:   PetscFree(pc->data);
1217:   return(0);
1218: }

1222: static PetscErrorCode PCApply_BJacobi_Multiproc(PC pc,Vec x,Vec y)
1223: {
1224:   PC_BJacobi           *jac = (PC_BJacobi*)pc->data;
1225:   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1226:   PetscErrorCode       ierr;
1227:   PetscScalar          *xarray,*yarray;

1230:   /* place x's and y's local arrays into xsub and ysub */
1231:   VecGetArray(x,&xarray);
1232:   VecGetArray(y,&yarray);
1233:   VecPlaceArray(mpjac->xsub,xarray);
1234:   VecPlaceArray(mpjac->ysub,yarray);

1236:   /* apply preconditioner on each matrix block */
1237:   KSPSolve(mpjac->ksp,mpjac->xsub,mpjac->ysub);

1239:   VecResetArray(mpjac->xsub);
1240:   VecResetArray(mpjac->ysub);
1241:   VecRestoreArray(x,&xarray);
1242:   VecRestoreArray(y,&yarray);
1243:   return(0);
1244: }

1249: static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC pc)
1250: {
1251:   PC_BJacobi           *jac = (PC_BJacobi*)pc->data;
1252:   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1253:   PetscErrorCode       ierr;
1254:   PetscInt             m,n;
1255:   MPI_Comm             comm = ((PetscObject)pc)->comm,subcomm=0;
1256:   const char           *prefix;

1259:   if (jac->n_local > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only a single block in a subcommunicator is supported");
1260:   jac->n_local = 1; /* currently only a single block is supported for a subcommunicator */
1261:   if (!pc->setupcalled) {
1262:     PetscNewLog(pc,PC_BJacobi_Multiproc,&mpjac);
1263:     jac->data = (void*)mpjac;

1265:     /* initialize datastructure mpjac */
1266:     if (!jac->psubcomm) {
1267:       /* Create default contiguous subcommunicatiors if user does not provide them */
1268:       PetscSubcommCreate(comm,&jac->psubcomm);
1269:       PetscSubcommSetNumber(jac->psubcomm,jac->n);
1270:       PetscSubcommSetType(jac->psubcomm,PETSC_SUBCOMM_CONTIGUOUS);
1271:       PetscLogObjectMemory(pc,sizeof(PetscSubcomm));
1272:     }
1273:     mpjac->psubcomm = jac->psubcomm;
1274:     subcomm         = mpjac->psubcomm->comm;

1276:     /* Get matrix blocks of pmat */
1277:     MatGetMultiProcBlock_MPIAIJ(pc->pmat,subcomm,&mpjac->submats);

1279:     /* create a new PC that processors in each subcomm have copy of */
1280:     KSPCreate(subcomm,&mpjac->ksp);
1281:     PetscObjectIncrementTabLevel((PetscObject)mpjac->ksp,(PetscObject)pc,1);
1282:     PetscLogObjectParent(pc,mpjac->ksp);
1283:     KSPSetOperators(mpjac->ksp,mpjac->submats,mpjac->submats,pc->flag);
1284:     KSPGetPC(mpjac->ksp,&mpjac->pc);

1286:     PCGetOptionsPrefix(pc,&prefix);
1287:     KSPSetOptionsPrefix(mpjac->ksp,prefix);
1288:     KSPAppendOptionsPrefix(mpjac->ksp,"sub_");
1289:     /*
1290:       PetscMPIInt rank,subsize,subrank;
1291:       MPI_Comm_rank(comm,&rank);
1292:       MPI_Comm_size(subcomm,&subsize);
1293:       MPI_Comm_rank(subcomm,&subrank);

1295:       MatGetLocalSize(mpjac->submats,&m,PETSC_NULL);
1296:       MatGetSize(mpjac->submats,&n,PETSC_NULL);
1297:       PetscSynchronizedPrintf(comm,"[%d], sub-size %d,sub-rank %d\n",rank,subsize,subrank);
1298:       PetscSynchronizedFlush(comm);
1299:     */

1301:     /* create dummy vectors xsub and ysub */
1302:     MatGetLocalSize(mpjac->submats,&m,&n);
1303:     VecCreateMPIWithArray(subcomm,n,PETSC_DECIDE,PETSC_NULL,&mpjac->xsub);
1304:     VecCreateMPIWithArray(subcomm,m,PETSC_DECIDE,PETSC_NULL,&mpjac->ysub);
1305:     PetscLogObjectParent(pc,mpjac->xsub);
1306:     PetscLogObjectParent(pc,mpjac->ysub);

1308:     pc->ops->reset   = PCReset_BJacobi_Multiproc;
1309:     pc->ops->destroy = PCDestroy_BJacobi_Multiproc;
1310:     pc->ops->apply   = PCApply_BJacobi_Multiproc;
1311:   }

1313:   if (pc->setupcalled && pc->flag == DIFFERENT_NONZERO_PATTERN) {
1314:     /* destroy old matrix blocks, then get new matrix blocks */
1315:     if (mpjac->submats) {
1316:       MatDestroy(&mpjac->submats);
1317:       subcomm = mpjac->psubcomm->comm;
1318:       MatGetMultiProcBlock_MPIAIJ(pc->pmat,subcomm,&mpjac->submats);
1319:       KSPSetOperators(mpjac->ksp,mpjac->submats,mpjac->submats,DIFFERENT_NONZERO_PATTERN);
1320:     }
1321:   }

1323:   if (pc->setfromoptionscalled){
1324:     KSPSetFromOptions(mpjac->ksp);
1325:   }
1326:   return(0);
1327: }