Actual source code: bjacobi.c

  1: #define PETSCKSP_DLL

  3: /*
  4:    Defines a block Jacobi preconditioner.
  5: */
 6:  #include include/private/matimpl.h
 7:  #include private/pcimpl.h
 8:  #include src/ksp/pc/impls/bjacobi/bjacobi.h

 10: static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC,Mat,Mat);
 11: static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC,Mat,Mat);

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

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

 31:   /* ----------
 32:       Determines the number of blocks assigned to each processor 
 33:   */

 35:   /*   local block count  given */
 36:   if (jac->n_local > 0 && jac->n < 0) {
 37:     MPI_Allreduce(&jac->n_local,&jac->n,1,MPIU_INT,MPI_SUM,pc->comm);
 38:     if (jac->l_lens) { /* check that user set these correctly */
 39:       sum = 0;
 40:       for (i=0; i<jac->n_local; i++) {
 41:         if (jac->l_lens[i]/bs*bs !=jac->l_lens[i]) {
 42:           SETERRQ(PETSC_ERR_ARG_SIZ,"Mat blocksize doesn't match block Jacobi layout");
 43:         }
 44:         sum += jac->l_lens[i];
 45:       }
 46:       if (sum != M) SETERRQ(PETSC_ERR_ARG_SIZ,"Local lens sent incorrectly");
 47:     } else {
 48:       PetscMalloc(jac->n_local*sizeof(PetscInt),&jac->l_lens);
 49:       for (i=0; i<jac->n_local; i++) {
 50:         jac->l_lens[i] = bs*((M/bs)/jac->n_local + (((M/bs) % jac->n_local) > i));
 51:       }
 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_ERR_ARG_SIZ,"Zero block not allowed");
 59:         if (jac->g_lens[i]/bs*bs != jac->g_lens[i]) {
 60:           SETERRQ(PETSC_ERR_ARG_SIZ,"Mat blocksize doesn't match block Jacobi layout");
 61:         }
 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_ERR_ARG_SIZ,"Global lens sent 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_ERR_ARG_SIZ,"Block sizes\n\
 80:                    used in PCBJacobiSetTotalBlocks()\n\
 81:                    are not compatible with parallel matrix layout");
 82:  start_1:
 83:         for (i=i_start; i<jac->n+1; i++) {
 84:           if (sum == end) { i_end = i; goto end_1; }
 85:           if (i < jac->n) sum += jac->g_lens[i];
 86:         }
 87:         SETERRQ(PETSC_ERR_ARG_SIZ,"Block sizes\n\
 88:                       used in PCBJacobiSetTotalBlocks()\n\
 89:                       are not compatible with parallel matrix layout");
 90:  end_1:
 91:         jac->n_local = i_end - i_start;
 92:         PetscMalloc(jac->n_local*sizeof(PetscInt),&jac->l_lens);
 93:         PetscMemcpy(jac->l_lens,jac->g_lens+i_start,jac->n_local*sizeof(PetscInt));
 94:       }
 95:     } else { /* no global blocks given, determine then using default layout */
 96:       jac->n_local = jac->n/size + ((jac->n % size) > rank);
 97:       PetscMalloc(jac->n_local*sizeof(PetscInt),&jac->l_lens);
 98:       for (i=0; i<jac->n_local; i++) {
 99:         jac->l_lens[i] = ((M/bs)/jac->n_local + (((M/bs) % jac->n_local) > i))*bs;
100:         if (!jac->l_lens[i]) SETERRQ(PETSC_ERR_ARG_SIZ,"Too many blocks given");
101:       }
102:     }
103:   } else if (jac->n < 0 && jac->n_local < 0) { /* no blocks given */
104:     jac->n         = size;
105:     jac->n_local   = 1;
106:     PetscMalloc(sizeof(PetscInt),&jac->l_lens);
107:     jac->l_lens[0] = M;
108:   }

110:   MPI_Comm_size(pc->comm,&size);
111:   PetscObjectQueryFunction((PetscObject)pc->mat,"MatGetDiagonalBlock_C",(void (**)(void))&f);
112:   if (size == 1 && !f) {
113:     mat  = pc->mat;
114:     pmat = pc->pmat;
115:   } else {
116:     PetscTruth iscopy;
117:     MatReuse   scall;

119:     if (jac->use_true_local) {
120:       scall = MAT_INITIAL_MATRIX;
121:       if (pc->setupcalled) {
122:         if (pc->flag == SAME_NONZERO_PATTERN) {
123:           if (jac->tp_mat) {
124:             scall = MAT_REUSE_MATRIX;
125:             mat   = jac->tp_mat;
126:           }
127:         } else {
128:           if (jac->tp_mat)  {
129:             MatDestroy(jac->tp_mat);
130:           }
131:         }
132:       }
133:       if (!f) {
134:         SETERRQ(PETSC_ERR_SUP,"This matrix does not support getting diagonal block");
135:       }
136:       (*f)(pc->mat,&iscopy,scall,&mat);
137:       /* make submatrix have same prefix as entire matrix */
138:       PetscObjectGetOptionsPrefix((PetscObject)pc->mat,&mprefix);
139:       PetscObjectSetOptionsPrefix((PetscObject)mat,mprefix);
140:       if (iscopy) {
141:         jac->tp_mat = mat;
142:       }
143:     }
144:     if (pc->pmat != pc->mat || !jac->use_true_local) {
145:       scall = MAT_INITIAL_MATRIX;
146:       if (pc->setupcalled) {
147:         if (pc->flag == SAME_NONZERO_PATTERN) {
148:           if (jac->tp_pmat) {
149:             scall = MAT_REUSE_MATRIX;
150:             pmat   = jac->tp_pmat;
151:           }
152:         } else {
153:           if (jac->tp_pmat)  {
154:             MatDestroy(jac->tp_pmat);
155:           }
156:         }
157:       }
158:       PetscObjectQueryFunction((PetscObject)pc->pmat,"MatGetDiagonalBlock_C",(void (**)(void))&f);
159:       if (!f) {
160:         const char *type;
161:         PetscObjectGetType((PetscObject) pc->pmat,&type);
162:         SETERRQ1(PETSC_ERR_SUP,"This matrix type, %s, does not support getting diagonal block", type);
163:       }
164:       (*f)(pc->pmat,&iscopy,scall,&pmat);
165:       /* make submatrix have same prefix as entire matrix */
166:       PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);
167:       PetscObjectSetOptionsPrefix((PetscObject)pmat,pprefix);
168:       if (iscopy) {
169:         jac->tp_pmat = pmat;
170:       }
171:     } else {
172:       pmat = mat;
173:     }
174:   }

176:   /* ------
177:      Setup code depends on the number of blocks 
178:   */
179:   if (jac->n_local == 1) {
180:     PCSetUp_BJacobi_Singleblock(pc,mat,pmat);
181:   } else {
182:     PCSetUp_BJacobi_Multiblock(pc,mat,pmat);
183:   }
184:   return(0);
185: }

187: /* Default destroy, if it has never been setup */
190: static PetscErrorCode PCDestroy_BJacobi(PC pc)
191: {
192:   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;

196:   PetscFree(jac->g_lens);
197:   PetscFree(jac->l_lens);
198:   PetscFree(jac);
199:   return(0);
200: }


205: static PetscErrorCode PCSetFromOptions_BJacobi(PC pc)
206: {
207:   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;
209:   PetscInt       blocks;
210:   PetscTruth     flg;

213:   PetscOptionsHead("Block Jacobi options");
214:     PetscOptionsInt("-pc_bjacobi_blocks","Total number of blocks","PCBJacobiSetTotalBlocks",jac->n,&blocks,&flg);
215:     if (flg) {
216:       PCBJacobiSetTotalBlocks(pc,blocks,PETSC_NULL);
217:     }
218:     PetscOptionsName("-pc_bjacobi_truelocal","Use the true matrix, not preconditioner matrix to define matrix vector product in sub-problems","PCBJacobiSetUseTrueLocal",&flg);
219:     if (flg) {
220:       PCBJacobiSetUseTrueLocal(pc);
221:     }
222:   PetscOptionsTail();
223:   return(0);
224: }

228: static PetscErrorCode PCView_BJacobi(PC pc,PetscViewer viewer)
229: {
230:   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;
232:   PetscMPIInt    rank;
233:   PetscInt       i;
234:   PetscTruth     iascii,isstring;
235:   PetscViewer    sviewer;

238:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
239:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
240:   if (iascii) {
241:     if (jac->use_true_local) {
242:       PetscViewerASCIIPrintf(viewer,"  block Jacobi: using true local matrix, number of blocks = %D\n",jac->n);
243:     }
244:     PetscViewerASCIIPrintf(viewer,"  block Jacobi: number of blocks = %D\n",jac->n);
245:     MPI_Comm_rank(pc->comm,&rank);
246:     if (jac->same_local_solves) {
247:       PetscViewerASCIIPrintf(viewer,"  Local solve is same for all blocks, in the following KSP and PC objects:\n");
248:       PetscViewerGetSingleton(viewer,&sviewer);
249:       if (!rank && jac->ksp) {
250:         PetscViewerASCIIPushTab(viewer);
251:         KSPView(jac->ksp[0],sviewer);
252:         PetscViewerASCIIPopTab(viewer);
253:       }
254:       PetscViewerRestoreSingleton(viewer,&sviewer);
255:     } else {
256:       PetscInt n_global;
257:       MPI_Allreduce(&jac->n_local,&n_global,1,MPIU_INT,MPI_MAX,pc->comm);
258:       PetscViewerASCIIPrintf(viewer,"  Local solve info for each block is in the following KSP and PC objects:\n");
259:       PetscViewerASCIISynchronizedPrintf(viewer,"[%d] number of local blocks = %D, first local block number = %D\n",
260:                    rank,jac->n_local,jac->first_local);
261:       PetscViewerASCIIPushTab(viewer);
262:       for (i=0; i<n_global; i++) {
263:         PetscViewerGetSingleton(viewer,&sviewer);
264:         if (i < jac->n_local) {
265:           PetscViewerASCIISynchronizedPrintf(viewer,"[%d] local block number %D\n",rank,i);
266:           KSPView(jac->ksp[i],sviewer);
267:           PetscViewerASCIISynchronizedPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");
268:         }
269:         PetscViewerRestoreSingleton(viewer,&sviewer);
270:       }
271:       PetscViewerASCIIPopTab(viewer);
272:       PetscViewerFlush(viewer);
273:     }
274:   } else if (isstring) {
275:     PetscViewerStringSPrintf(viewer," blks=%D",jac->n);
276:     PetscViewerGetSingleton(viewer,&sviewer);
277:     if (jac->ksp) {KSPView(jac->ksp[0],sviewer);}
278:     PetscViewerRestoreSingleton(viewer,&sviewer);
279:   } else {
280:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for block Jacobi",((PetscObject)viewer)->type_name);
281:   }
282:   return(0);
283: }

285: /* -------------------------------------------------------------------------------------*/

290: PetscErrorCode  PCBJacobiSetUseTrueLocal_BJacobi(PC pc)
291: {
292:   PC_BJacobi   *jac;

295:   jac                 = (PC_BJacobi*)pc->data;
296:   jac->use_true_local = PETSC_TRUE;
297:   return(0);
298: }

304: PetscErrorCode  PCBJacobiGetSubKSP_BJacobi(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp)
305: {
306:   PC_BJacobi   *jac = (PC_BJacobi*)pc->data;;

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

311:   if (n_local)     *n_local     = jac->n_local;
312:   if (first_local) *first_local = jac->first_local;
313:   *ksp                          = jac->ksp;
314:   jac->same_local_solves        = PETSC_FALSE; /* Assume that local solves are now different;
315:                                                   not necessarily true though!  This flag is 
316:                                                   used only for PCView_BJacobi() */
317:   return(0);
318: }

324: PetscErrorCode  PCBJacobiSetTotalBlocks_BJacobi(PC pc,PetscInt blocks,PetscInt *lens)
325: {
326:   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;


331:   if (pc->setupcalled > 0 && jac->n!=blocks) SETERRQ(PETSC_ERR_ORDER,"Cannot alter number of blocks after PCSetUp()/KSPSetUp() has been called");
332:   jac->n = blocks;
333:   if (!lens) {
334:     jac->g_lens = 0;
335:   } else {
336:     PetscMalloc(blocks*sizeof(PetscInt),&jac->g_lens);
337:     PetscLogObjectMemory(pc,blocks*sizeof(PetscInt));
338:     PetscMemcpy(jac->g_lens,lens,blocks*sizeof(PetscInt));
339:   }
340:   return(0);
341: }

347: PetscErrorCode  PCBJacobiGetTotalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[])
348: {
349:   PC_BJacobi *jac = (PC_BJacobi*) pc->data;

352:   *blocks = jac->n;
353:   if (lens) *lens = jac->g_lens;
354:   return(0);
355: }

361: PetscErrorCode  PCBJacobiSetLocalBlocks_BJacobi(PC pc,PetscInt blocks,const PetscInt lens[])
362: {
363:   PC_BJacobi     *jac;

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

369:   jac->n_local = blocks;
370:   if (!lens) {
371:     jac->l_lens = 0;
372:   } else {
373:     PetscMalloc(blocks*sizeof(PetscInt),&jac->l_lens);
374:     PetscLogObjectMemory(pc,blocks*sizeof(PetscInt));
375:     PetscMemcpy(jac->l_lens,lens,blocks*sizeof(PetscInt));
376:   }
377:   return(0);
378: }

384: PetscErrorCode  PCBJacobiGetLocalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[])
385: {
386:   PC_BJacobi *jac = (PC_BJacobi*) pc->data;

389:   *blocks = jac->n_local;
390:   if (lens) *lens = jac->l_lens;
391:   return(0);
392: }

395: /* -------------------------------------------------------------------------------------*/

399: /*@
400:    PCBJacobiSetUseTrueLocal - Sets a flag to indicate that the block 
401:    problem is associated with the linear system matrix instead of the
402:    default (where it is associated with the preconditioning matrix).
403:    That is, if the local system is solved iteratively then it iterates
404:    on the block from the matrix using the block from the preconditioner
405:    as the preconditioner for the local block.

407:    Collective on PC

409:    Input Parameters:
410: .  pc - the preconditioner context

412:    Options Database Key:
413: .  -pc_bjacobi_truelocal - Activates PCBJacobiSetUseTrueLocal()

415:    Notes:
416:    For the common case in which the preconditioning and linear 
417:    system matrices are identical, this routine is unnecessary.

419:    Level: intermediate

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

423: .seealso: PCSetOperators(), PCBJacobiSetLocalBlocks()
424: @*/
425: PetscErrorCode  PCBJacobiSetUseTrueLocal(PC pc)
426: {
427:   PetscErrorCode ierr,(*f)(PC);

431:   PetscObjectQueryFunction((PetscObject)pc,"PCBJacobiSetUseTrueLocal_C",(void (**)(void))&f);
432:   if (f) {
433:     (*f)(pc);
434:   }

436:   return(0);
437: }

441: /*@C
442:    PCBJacobiGetSubKSP - Gets the local KSP contexts for all blocks on
443:    this processor.
444:    
445:    Note Collective

447:    Input Parameter:
448: .  pc - the preconditioner context

450:    Output Parameters:
451: +  n_local - the number of blocks on this processor, or PETSC_NULL
452: .  first_local - the global number of the first block on this processor, or PETSC_NULL
453: -  ksp - the array of KSP contexts

455:    Notes:  
456:    After PCBJacobiGetSubKSP() the array of KSP contexts is not to be freed.
457:    
458:    Currently for some matrix implementations only 1 block per processor 
459:    is supported.
460:    
461:    You must call KSPSetUp() or PCSetUp() before calling PCBJacobiGetSubKSP().

463:    Level: advanced

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

467: .seealso: PCBJacobiGetSubKSP()
468: @*/
469: PetscErrorCode  PCBJacobiGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
470: {
471:   PetscErrorCode ierr,(*f)(PC,PetscInt *,PetscInt *,KSP **);

475:   PetscObjectQueryFunction((PetscObject)pc,"PCBJacobiGetSubKSP_C",(void (**)(void))&f);
476:   if (f) {
477:     (*f)(pc,n_local,first_local,ksp);
478:   } else {
479:     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subsolvers for this preconditioner");
480:   }
481:   return(0);
482: }

486: /*@
487:    PCBJacobiSetTotalBlocks - Sets the global number of blocks for the block
488:    Jacobi preconditioner.

490:    Collective on PC

492:    Input Parameters:
493: +  pc - the preconditioner context
494: .  blocks - the number of blocks
495: -  lens - [optional] integer array containing the size of each block

497:    Options Database Key:
498: .  -pc_bjacobi_blocks <blocks> - Sets the number of global blocks

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

504:    Level: intermediate

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

508: .seealso: PCBJacobiSetUseTrueLocal(), PCBJacobiSetLocalBlocks()
509: @*/
510: PetscErrorCode  PCBJacobiSetTotalBlocks(PC pc,PetscInt blocks,const PetscInt lens[])
511: {
512:   PetscErrorCode ierr,(*f)(PC,PetscInt,const PetscInt[]);

516:   if (blocks <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Must have positive blocks");
517:   PetscObjectQueryFunction((PetscObject)pc,"PCBJacobiSetTotalBlocks_C",(void (**)(void))&f);
518:   if (f) {
519:     (*f)(pc,blocks,lens);
520:   }
521:   return(0);
522: }

526: /*@C
527:    PCBJacobiGetTotalBlocks - Gets the global number of blocks for the block
528:    Jacobi preconditioner.

530:    Collective on PC

532:    Input Parameter:
533: .  pc - the preconditioner context

535:    Output parameters:
536: +  blocks - the number of blocks
537: -  lens - integer array containing the size of each block

539:    Level: intermediate

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

543: .seealso: PCBJacobiSetUseTrueLocal(), PCBJacobiGetLocalBlocks()
544: @*/
545: PetscErrorCode  PCBJacobiGetTotalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
546: {
547:   PetscErrorCode ierr,(*f)(PC,PetscInt*, const PetscInt *[]);

552:   PetscObjectQueryFunction((PetscObject)pc,"PCBJacobiGetTotalBlocks_C",(void (**)(void))&f);
553:   if (f) {
554:     (*f)(pc,blocks,lens);
555:   }
556:   return(0);
557: }
558: 
561: /*@
562:    PCBJacobiSetLocalBlocks - Sets the local number of blocks for the block
563:    Jacobi preconditioner.

565:    Not Collective

567:    Input Parameters:
568: +  pc - the preconditioner context
569: .  blocks - the number of blocks
570: -  lens - [optional] integer array containing size of each block

572:    Note:  
573:    Currently only a limited number of blocking configurations are supported.

575:    Level: intermediate

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

579: .seealso: PCBJacobiSetUseTrueLocal(), PCBJacobiSetTotalBlocks()
580: @*/
581: PetscErrorCode  PCBJacobiSetLocalBlocks(PC pc,PetscInt blocks,const PetscInt lens[])
582: {
583:   PetscErrorCode ierr,(*f)(PC,PetscInt,const PetscInt []);

587:   if (blocks < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Must have nonegative blocks");
588:   PetscObjectQueryFunction((PetscObject)pc,"PCBJacobiSetLocalBlocks_C",(void (**)(void))&f);
589:   if (f) {
590:     (*f)(pc,blocks,lens);
591:   }
592:   return(0);
593: }
594: 
597: /*@C
598:    PCBJacobiGetLocalBlocks - Gets the local number of blocks for the block
599:    Jacobi preconditioner.

601:    Not Collective

603:    Input Parameters:
604: +  pc - the preconditioner context
605: .  blocks - the number of blocks
606: -  lens - [optional] integer array containing size of each block

608:    Note:  
609:    Currently only a limited number of blocking configurations are supported.

611:    Level: intermediate

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

615: .seealso: PCBJacobiSetUseTrueLocal(), PCBJacobiGetTotalBlocks()
616: @*/
617: PetscErrorCode  PCBJacobiGetLocalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
618: {
619:   PetscErrorCode ierr,(*f)(PC,PetscInt*, const PetscInt *[]);

624:   PetscObjectQueryFunction((PetscObject)pc,"PCBJacobiGetLocalBlocks_C",(void (**)(void))&f);
625:   if (f) {
626:     (*f)(pc,blocks,lens);
627:   }
628:   return(0);
629: }

631: /* -----------------------------------------------------------------------------------*/

633: /*MC
634:    PCBJACOBI - Use block Jacobi preconditioning, each block is (approximately) solved with 
635:            its own KSP object.

637:    Options Database Keys:
638: .  -pc_bjacobi_truelocal - Activates PCBJacobiSetUseTrueLocal()

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

643:      To set options on the solvers for each block append -sub_ to all the KSP, KSP, and PC
644:         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
645:         
646:      To set the options on the solvers separate for each block call PCBJacobiGetSubKSP()
647:          and set the options directly on the resulting KSP object (you can access its PC
648:          KSPGetPC())

650:    Level: beginner

652:    Concepts: block Jacobi

654: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
655:            PCASM, PCBJacobiSetUseTrueLocal(), PCBJacobiGetSubKSP(), PCBJacobiSetTotalBlocks(),
656:            PCBJacobiSetLocalBlocks(), PCSetModifySubmatrices()
657: M*/

662: PetscErrorCode  PCCreate_BJacobi(PC pc)
663: {
665:   PetscMPIInt    rank;
666:   PC_BJacobi     *jac;

669:   PetscNew(PC_BJacobi,&jac);
670:   PetscLogObjectMemory(pc,sizeof(PC_BJacobi));
671:   MPI_Comm_rank(pc->comm,&rank);
672:   pc->ops->apply              = 0;
673:   pc->ops->applytranspose     = 0;
674:   pc->ops->setup              = PCSetUp_BJacobi;
675:   pc->ops->destroy            = PCDestroy_BJacobi;
676:   pc->ops->setfromoptions     = PCSetFromOptions_BJacobi;
677:   pc->ops->view               = PCView_BJacobi;
678:   pc->ops->applyrichardson    = 0;

680:   pc->data               = (void*)jac;
681:   jac->n                 = -1;
682:   jac->n_local           = -1;
683:   jac->first_local       = rank;
684:   jac->ksp              = 0;
685:   jac->use_true_local    = PETSC_FALSE;
686:   jac->same_local_solves = PETSC_TRUE;
687:   jac->g_lens            = 0;
688:   jac->l_lens            = 0;
689:   jac->tp_mat            = 0;
690:   jac->tp_pmat           = 0;

692:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiSetUseTrueLocal_C",
693:                     "PCBJacobiSetUseTrueLocal_BJacobi",
694:                     PCBJacobiSetUseTrueLocal_BJacobi);
695:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiGetSubKSP_C","PCBJacobiGetSubKSP_BJacobi",
696:                     PCBJacobiGetSubKSP_BJacobi);
697:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiSetTotalBlocks_C","PCBJacobiSetTotalBlocks_BJacobi",
698:                     PCBJacobiSetTotalBlocks_BJacobi);
699:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiGetTotalBlocks_C","PCBJacobiGetTotalBlocks_BJacobi",
700:                     PCBJacobiGetTotalBlocks_BJacobi);
701:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiSetLocalBlocks_C","PCBJacobiSetLocalBlocks_BJacobi",
702:                     PCBJacobiSetLocalBlocks_BJacobi);
703:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiGetLocalBlocks_C","PCBJacobiGetLocalBlocks_BJacobi",
704:                     PCBJacobiGetLocalBlocks_BJacobi);

706:   return(0);
707: }

710: /* --------------------------------------------------------------------------------------------*/
711: /*
712:         These are for a single block per processor; works for AIJ, BAIJ; Seq and MPI
713: */
716: PetscErrorCode PCDestroy_BJacobi_Singleblock(PC pc)
717: {
718:   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
719:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
720:   PetscErrorCode         ierr;

723:   /*
724:         If the on processor block had to be generated via a MatGetDiagonalBlock()
725:      that creates a copy (for example MPIBDiag matrices do), this frees the space
726:   */
727:   if (jac->tp_mat) {
728:     MatDestroy(jac->tp_mat);
729:   }
730:   if (jac->tp_pmat) {
731:     MatDestroy(jac->tp_pmat);
732:   }

734:   KSPDestroy(jac->ksp[0]);
735:   PetscFree(jac->ksp);
736:   VecDestroy(bjac->x);
737:   VecDestroy(bjac->y);
738:   PetscFree(jac->l_lens);
739:   PetscFree(jac->g_lens);
740:   PetscFree(bjac);
741:   PetscFree(jac);
742:   return(0);
743: }

747: PetscErrorCode PCSetUpOnBlocks_BJacobi_Singleblock(PC pc)
748: {
750:   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;

753:   KSPSetUp(jac->ksp[0]);
754:   return(0);
755: }

759: PetscErrorCode PCApply_BJacobi_Singleblock(PC pc,Vec x,Vec y)
760: {
761:   PetscErrorCode         ierr;
762:   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
763:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
764:   PetscScalar            *x_array,*y_array;

767:   /* 
768:       The VecPlaceArray() is to avoid having to copy the 
769:     y vector into the bjac->x vector. The reason for 
770:     the bjac->x vector is that we need a sequential vector
771:     for the sequential solve.
772:   */
773:   VecGetArray(x,&x_array);
774:   VecGetArray(y,&y_array);
775:   VecPlaceArray(bjac->x,x_array);
776:   VecPlaceArray(bjac->y,y_array);
777:   KSPSolve(jac->ksp[0],bjac->x,bjac->y);
778:   VecResetArray(bjac->x);
779:   VecResetArray(bjac->y);
780:   VecRestoreArray(x,&x_array);
781:   VecRestoreArray(y,&y_array);
782:   return(0);
783: }

787: PetscErrorCode PCApplySymmetricLeft_BJacobi_Singleblock(PC pc,Vec x,Vec y)
788: {
789:   PetscErrorCode         ierr;
790:   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
791:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
792:   PetscScalar            *x_array,*y_array;
793:   PC                     subpc;

796:   /* 
797:       The VecPlaceArray() is to avoid having to copy the 
798:     y vector into the bjac->x vector. The reason for 
799:     the bjac->x vector is that we need a sequential vector
800:     for the sequential solve.
801:   */
802:   VecGetArray(x,&x_array);
803:   VecGetArray(y,&y_array);
804:   VecPlaceArray(bjac->x,x_array);
805:   VecPlaceArray(bjac->y,y_array);

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

810:   KSPGetPC(jac->ksp[0],&subpc);
811:   PCApplySymmetricLeft(subpc,bjac->x,bjac->y);
812:   VecResetArray(bjac->x);
813:   VecResetArray(bjac->y);

815:   VecRestoreArray(x,&x_array);
816:   VecRestoreArray(y,&y_array);
817:   return(0);
818: }

822: PetscErrorCode PCApplySymmetricRight_BJacobi_Singleblock(PC pc,Vec x,Vec y)
823: {
824:   PetscErrorCode         ierr;
825:   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
826:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
827:   PetscScalar            *x_array,*y_array;
828:   PC                     subpc;

831:   /* 
832:       The VecPlaceArray() is to avoid having to copy the 
833:     y vector into the bjac->x vector. The reason for 
834:     the bjac->x vector is that we need a sequential vector
835:     for the sequential solve.
836:   */
837:   VecGetArray(x,&x_array);
838:   VecGetArray(y,&y_array);
839:   VecPlaceArray(bjac->x,x_array);
840:   VecPlaceArray(bjac->y,y_array);

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

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

848:   VecRestoreArray(x,&x_array);
849:   VecRestoreArray(y,&y_array);
850:   return(0);
851: }

855: PetscErrorCode PCApplyTranspose_BJacobi_Singleblock(PC pc,Vec x,Vec y)
856: {
857:   PetscErrorCode         ierr;
858:   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
859:   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
860:   PetscScalar            *x_array,*y_array;

863:   /* 
864:       The VecPlaceArray() is to avoid having to copy the 
865:     y vector into the bjac->x vector. The reason for 
866:     the bjac->x vector is that we need a sequential vector
867:     for the sequential solve.
868:   */
869:   VecGetArray(x,&x_array);
870:   VecGetArray(y,&y_array);
871:   VecPlaceArray(bjac->x,x_array);
872:   VecPlaceArray(bjac->y,y_array);
873:   KSPSolveTranspose(jac->ksp[0],bjac->x,bjac->y);
874:   VecResetArray(bjac->x);
875:   VecResetArray(bjac->y);
876:   VecRestoreArray(x,&x_array);
877:   VecRestoreArray(y,&y_array);
878:   return(0);
879: }

883: static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC pc,Mat mat,Mat pmat)
884: {
885:   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
886:   PetscErrorCode         ierr;
887:   PetscInt               m;
888:   KSP                    ksp;
889:   Vec                    x,y;
890:   PC_BJacobi_Singleblock *bjac;
891:   PC                     subpc;
892:   PetscTruth             wasSetup;


896:   /* set default direct solver with no Krylov method */
897:   if (!pc->setupcalled) {
898:     const char *prefix;
899:     wasSetup = PETSC_FALSE;
900:     KSPCreate(PETSC_COMM_SELF,&ksp);
901:     PetscLogObjectParent(pc,ksp);
902:     KSPSetType(ksp,KSPPREONLY);
903:     KSPGetPC(ksp,&subpc);
904:     PCGetOptionsPrefix(pc,&prefix);
905:     KSPSetOptionsPrefix(ksp,prefix);
906:     KSPAppendOptionsPrefix(ksp,"sub_");
907:     /*
908:       The reason we need to generate these vectors is to serve 
909:       as the right-hand side and solution vector for the solve on the 
910:       block. We do not need to allocate space for the vectors since
911:       that is provided via VecPlaceArray() just before the call to 
912:       KSPSolve() on the block.
913:     */
914:     MatGetSize(pmat,&m,&m);
915:     VecCreateSeqWithArray(PETSC_COMM_SELF,m,PETSC_NULL,&x);
916:     VecCreateSeqWithArray(PETSC_COMM_SELF,m,PETSC_NULL,&y);
917:     PetscLogObjectParent(pc,x);
918:     PetscLogObjectParent(pc,y);

920:     pc->ops->destroy             = PCDestroy_BJacobi_Singleblock;
921:     pc->ops->apply               = PCApply_BJacobi_Singleblock;
922:     pc->ops->applysymmetricleft  = PCApplySymmetricLeft_BJacobi_Singleblock;
923:     pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Singleblock;
924:     pc->ops->applytranspose      = PCApplyTranspose_BJacobi_Singleblock;
925:     pc->ops->setuponblocks       = PCSetUpOnBlocks_BJacobi_Singleblock;

927:     PetscMalloc(sizeof(PC_BJacobi_Singleblock),&bjac);
928:     PetscLogObjectMemory(pc,sizeof(PC_BJacobi_Singleblock));
929:     bjac->x      = x;
930:     bjac->y      = y;

932:     PetscMalloc(sizeof(KSP),&jac->ksp);
933:     jac->ksp[0] = ksp;
934:     jac->data    = (void*)bjac;
935:   } else {
936:     wasSetup = PETSC_TRUE;
937:     ksp = jac->ksp[0];
938:     bjac = (PC_BJacobi_Singleblock *)jac->data;
939:   }
940:   if (jac->use_true_local) {
941:     KSPSetOperators(ksp,mat,pmat,pc->flag);
942:   }  else {
943:     KSPSetOperators(ksp,pmat,pmat,pc->flag);
944:   }
945:   if (!wasSetup) {
946:     KSPSetFromOptions(ksp);
947:   }
948:   return(0);
949: }

951: /* ---------------------------------------------------------------------------------------------*/

955: PetscErrorCode PCDestroy_BJacobi_Multiblock(PC pc)
956: {
957:   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
958:   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
959:   PetscErrorCode        ierr;
960:   PetscInt              i;

963:   MatDestroyMatrices(jac->n_local,&bjac->pmat);
964:   if (jac->use_true_local) {
965:     MatDestroyMatrices(jac->n_local,&bjac->mat);
966:   }

968:   /*
969:         If the on processor block had to be generated via a MatGetDiagonalBlock()
970:      that creates a copy (for example MPIBDiag matrices do), this frees the space
971:   */
972:   if (jac->tp_mat) {
973:     MatDestroy(jac->tp_mat);
974:   }
975:   if (jac->tp_pmat) {
976:     MatDestroy(jac->tp_pmat);
977:   }

979:   for (i=0; i<jac->n_local; i++) {
980:     KSPDestroy(jac->ksp[i]);
981:     VecDestroy(bjac->x[i]);
982:     VecDestroy(bjac->y[i]);
983:     ISDestroy(bjac->is[i]);
984:   }
985:   PetscFree(jac->ksp);
986:   PetscFree2(bjac->x,bjac->y);
987:   PetscFree(bjac->starts);
988:   PetscFree(bjac->is);
989:   PetscFree(bjac);
990:   PetscFree(jac->l_lens);
991:   PetscFree(jac->g_lens);
992:   PetscFree(jac);
993:   return(0);
994: }

998: PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiblock(PC pc)
999: {
1000:   PC_BJacobi     *jac = (PC_BJacobi*)pc->data;
1002:   PetscInt       i,n_local = jac->n_local;

1005:   for (i=0; i<n_local; i++) {
1006:     KSPSetUp(jac->ksp[i]);
1007:   }
1008:   return(0);
1009: }

1011: /*
1012:       Preconditioner for block Jacobi 
1013: */
1016: PetscErrorCode PCApply_BJacobi_Multiblock(PC pc,Vec x,Vec y)
1017: {
1018:   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
1019:   PetscErrorCode        ierr;
1020:   PetscInt              i,n_local = jac->n_local;
1021:   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
1022:   PetscScalar           *xin,*yin;
1023:   static PetscTruth     flag = PETSC_TRUE;
1024: #if defined (PETSC_USE_LOG)
1025:   static PetscEvent     SUBKspSolve;
1026: #endif
1028:   if (flag) {
1030:     flag = PETSC_FALSE;
1031:   }
1032:   VecGetArray(x,&xin);
1033:   VecGetArray(y,&yin);
1034:   for (i=0; i<n_local; i++) {
1035:     /* 
1036:        To avoid copying the subvector from x into a workspace we instead 
1037:        make the workspace vector array point to the subpart of the array of
1038:        the global vector.
1039:     */
1040:     VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);
1041:     VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);

1044:     KSPSolve(jac->ksp[i],bjac->x[i],bjac->y[i]);

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

1055: /*
1056:       Preconditioner for block Jacobi 
1057: */
1060: PetscErrorCode PCApplyTranspose_BJacobi_Multiblock(PC pc,Vec x,Vec y)
1061: {
1062:   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
1063:   PetscErrorCode        ierr;
1064:   PetscInt              i,n_local = jac->n_local;
1065:   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
1066:   PetscScalar           *xin,*yin;
1067:   static PetscTruth     flag = PETSC_TRUE;
1068: #if defined (PETSC_USE_LOG)
1069:   static PetscEvent     SUBKspSolve;
1070: #endif

1073:   if (flag) {
1075:     flag = PETSC_FALSE;
1076:   }
1077:   VecGetArray(x,&xin);
1078:   VecGetArray(y,&yin);
1079:   for (i=0; i<n_local; i++) {
1080:     /* 
1081:        To avoid copying the subvector from x into a workspace we instead 
1082:        make the workspace vector array point to the subpart of the array of
1083:        the global vector.
1084:     */
1085:     VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);
1086:     VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);

1089:     KSPSolveTranspose(jac->ksp[i],bjac->x[i],bjac->y[i]);
1091:   }
1092:   VecRestoreArray(x,&xin);
1093:   VecRestoreArray(y,&yin);
1094:   return(0);
1095: }

1099: static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC pc,Mat mat,Mat pmat)
1100: {
1101:   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
1102:   PetscErrorCode         ierr;
1103:   PetscInt               m,n_local,N,M,start,i;
1104:   const char             *prefix,*pprefix,*mprefix;
1105:   KSP                    ksp;
1106:   Vec                    x,y;
1107:   PC_BJacobi_Multiblock  *bjac = (PC_BJacobi_Multiblock*)jac->data;
1108:   PC                     subpc;
1109:   IS                     is;
1110:   MatReuse               scall = MAT_REUSE_MATRIX;

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

1115:   n_local = jac->n_local;

1117:   if (jac->use_true_local) {
1118:     if (mat->type != pmat->type) SETERRQ(PETSC_ERR_ARG_INCOMP,"Matrices not of same type");
1119:   }

1121:   if (!pc->setupcalled) {
1122:     scall                  = MAT_INITIAL_MATRIX;
1123:     pc->ops->destroy       = PCDestroy_BJacobi_Multiblock;
1124:     pc->ops->apply         = PCApply_BJacobi_Multiblock;
1125:     pc->ops->applytranspose= PCApplyTranspose_BJacobi_Multiblock;
1126:     pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiblock;

1128:     PetscMalloc(sizeof(PC_BJacobi_Multiblock),&bjac);
1129:     PetscLogObjectMemory(pc,sizeof(PC_BJacobi_Multiblock));
1130:     PetscMalloc(n_local*sizeof(KSP),&jac->ksp);
1131:     PetscLogObjectMemory(pc,sizeof(n_local*sizeof(KSP)));
1132:     PetscMalloc2(n_local,Vec,&bjac->x,n_local,Vec,&bjac->y);
1133:     PetscMalloc(n_local*sizeof(PetscScalar),&bjac->starts);
1134:     PetscLogObjectMemory(pc,sizeof(n_local*sizeof(PetscScalar)));
1135: 
1136:     jac->data    = (void*)bjac;
1137:     PetscMalloc(n_local*sizeof(IS),&bjac->is);
1138:     PetscLogObjectMemory(pc,sizeof(n_local*sizeof(IS)));

1140:     start = 0;
1141:     for (i=0; i<n_local; i++) {
1142:       KSPCreate(PETSC_COMM_SELF,&ksp);
1143:       PetscLogObjectParent(pc,ksp);
1144:       KSPSetType(ksp,KSPPREONLY);
1145:       KSPGetPC(ksp,&subpc);
1146:       PCGetOptionsPrefix(pc,&prefix);
1147:       KSPSetOptionsPrefix(ksp,prefix);
1148:       KSPAppendOptionsPrefix(ksp,"sub_");

1150:       m = jac->l_lens[i];

1152:       /*
1153:       The reason we need to generate these vectors is to serve 
1154:       as the right-hand side and solution vector for the solve on the 
1155:       block. We do not need to allocate space for the vectors since
1156:       that is provided via VecPlaceArray() just before the call to 
1157:       KSPSolve() on the block.

1159:       */
1160:       VecCreateSeq(PETSC_COMM_SELF,m,&x);
1161:       VecCreateSeqWithArray(PETSC_COMM_SELF,m,PETSC_NULL,&y);
1162:       PetscLogObjectParent(pc,x);
1163:       PetscLogObjectParent(pc,y);
1164:       bjac->x[i]      = x;
1165:       bjac->y[i]      = y;
1166:       bjac->starts[i] = start;
1167:       jac->ksp[i]    = ksp;

1169:       ISCreateStride(PETSC_COMM_SELF,m,start,1,&is);
1170:       bjac->is[i] = is;
1171:       PetscLogObjectParent(pc,is);

1173:       start += m;
1174:     }
1175:   } else {
1176:     bjac = (PC_BJacobi_Multiblock*)jac->data;
1177:     /* 
1178:        Destroy the blocks from the previous iteration
1179:     */
1180:     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
1181:       MatDestroyMatrices(n_local,&bjac->pmat);
1182:       if (jac->use_true_local) {
1183:         MatDestroyMatrices(n_local,&bjac->mat);
1184:       }
1185:       scall = MAT_INITIAL_MATRIX;
1186:     }
1187:   }

1189:   MatGetSubMatrices(pmat,n_local,bjac->is,bjac->is,scall,&bjac->pmat);
1190:   if (jac->use_true_local) {
1191:     PetscObjectGetOptionsPrefix((PetscObject)mat,&mprefix);
1192:     MatGetSubMatrices(mat,n_local,bjac->is,bjac->is,scall,&bjac->mat);
1193:   }
1194:   /* Return control to the user so that the submatrices can be modified (e.g., to apply
1195:      different boundary conditions for the submatrices than for the global problem) */
1196:   PCModifySubMatrices(pc,n_local,bjac->is,bjac->is,bjac->pmat,pc->modifysubmatricesP);

1198:   PetscObjectGetOptionsPrefix((PetscObject)pmat,&pprefix);
1199:   for (i=0; i<n_local; i++) {
1200:     PetscLogObjectParent(pc,bjac->pmat[i]);
1201:     PetscObjectSetOptionsPrefix((PetscObject)bjac->pmat[i],pprefix);
1202:     if (jac->use_true_local) {
1203:       PetscLogObjectParent(pc,bjac->mat[i]);
1204:       PetscObjectSetOptionsPrefix((PetscObject)bjac->mat[i],mprefix);
1205:       KSPSetOperators(jac->ksp[i],bjac->mat[i],bjac->pmat[i],pc->flag);
1206:     } else {
1207:       KSPSetOperators(jac->ksp[i],bjac->pmat[i],bjac->pmat[i],pc->flag);
1208:     }
1209:     KSPSetFromOptions(jac->ksp[i]);
1210:   }

1212:   return(0);
1213: }