Actual source code: bjacobi.c

  1: #define PETSCKSP_DLL

  3: /*
  4:    Defines a block Jacobi preconditioner.
  5: */
 6:  #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(((PetscObject)pc)->comm,&rank);
 27:   MPI_Comm_size(((PetscObject)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,((PetscObject)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:   }
109:   if (jac->n_local < 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Number of blocks is less than number of processors");

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

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

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

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

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


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

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

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

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

287: /* -------------------------------------------------------------------------------------*/

292: PetscErrorCode  PCBJacobiSetUseTrueLocal_BJacobi(PC pc)
293: {
294:   PC_BJacobi   *jac;

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

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

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

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

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


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

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

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

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

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

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

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

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

397: /* -------------------------------------------------------------------------------------*/

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

409:    Collective on PC

411:    Input Parameters:
412: .  pc - the preconditioner context

414:    Options Database Key:
415: .  -pc_bjacobi_truelocal - Activates PCBJacobiSetUseTrueLocal()

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

421:    Level: intermediate

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

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

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

438:   return(0);
439: }

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

449:    Input Parameter:
450: .  pc - the preconditioner context

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

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

465:    Level: advanced

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

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

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

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

492:    Collective on PC

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

499:    Options Database Key:
500: .  -pc_bjacobi_blocks <blocks> - Sets the number of global blocks

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

506:    Level: intermediate

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

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

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

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

532:    Collective on PC

534:    Input Parameter:
535: .  pc - the preconditioner context

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

541:    Level: intermediate

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

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

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

567:    Not Collective

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

574:    Note:  
575:    Currently only a limited number of blocking configurations are supported.

577:    Level: intermediate

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

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

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

603:    Not Collective

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

610:    Note:  
611:    Currently only a limited number of blocking configurations are supported.

613:    Level: intermediate

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

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

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

633: /* -----------------------------------------------------------------------------------*/

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

639:    Options Database Keys:
640: .  -pc_bjacobi_truelocal - Activates PCBJacobiSetUseTrueLocal()

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

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

652:    Level: beginner

654:    Concepts: block Jacobi

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

664: PetscErrorCode  PCCreate_BJacobi(PC pc)
665: {
667:   PetscMPIInt    rank;
668:   PC_BJacobi     *jac;

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

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

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

707:   return(0);
708: }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

884: static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC pc,Mat mat,Mat pmat)
885: {
886:   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
887:   PetscErrorCode         ierr;
888:   PetscInt               m;
889:   KSP                    ksp;
890:   Vec                    x,y;
891:   PC_BJacobi_Singleblock *bjac;
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:     PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);
902:     PetscLogObjectParent(pc,ksp);
903:     KSPSetType(ksp,KSPPREONLY);
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:     PetscNewLog(pc,PC_BJacobi_Singleblock,&bjac);
928:     bjac->x      = x;
929:     bjac->y      = y;

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

950: /* ---------------------------------------------------------------------------------------------*/

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

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

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

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

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

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

1010: /*
1011:       Preconditioner for block Jacobi 
1012: */
1015: PetscErrorCode PCApply_BJacobi_Multiblock(PC pc,Vec x,Vec y)
1016: {
1017:   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
1018:   PetscErrorCode        ierr;
1019:   PetscInt              i,n_local = jac->n_local;
1020:   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
1021:   PetscScalar           *xin,*yin;

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

1035:     PetscLogEventBegin(PC_SetUpOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);
1036:     KSPSolve(jac->ksp[i],bjac->x[i],bjac->y[i]);
1037:     PetscLogEventEnd(PC_SetUpOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);

1039:     VecResetArray(bjac->x[i]);
1040:     VecResetArray(bjac->y[i]);
1041:   }
1042:   VecRestoreArray(x,&xin);
1043:   VecRestoreArray(y,&yin);
1044:   return(0);
1045: }

1047: /*
1048:       Preconditioner for block Jacobi 
1049: */
1052: PetscErrorCode PCApplyTranspose_BJacobi_Multiblock(PC pc,Vec x,Vec y)
1053: {
1054:   PC_BJacobi            *jac = (PC_BJacobi*)pc->data;
1055:   PetscErrorCode        ierr;
1056:   PetscInt              i,n_local = jac->n_local;
1057:   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
1058:   PetscScalar           *xin,*yin;

1061:   VecGetArray(x,&xin);
1062:   VecGetArray(y,&yin);
1063:   for (i=0; i<n_local; i++) {
1064:     /* 
1065:        To avoid copying the subvector from x into a workspace we instead 
1066:        make the workspace vector array point to the subpart of the array of
1067:        the global vector.
1068:     */
1069:     VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);
1070:     VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);

1072:     PetscLogEventBegin(PC_ApplyTransposeOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);
1073:     KSPSolveTranspose(jac->ksp[i],bjac->x[i],bjac->y[i]);
1074:     PetscLogEventEnd(PC_ApplyTransposeOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);
1075:   }
1076:   VecRestoreArray(x,&xin);
1077:   VecRestoreArray(y,&yin);
1078:   return(0);
1079: }

1083: static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC pc,Mat mat,Mat pmat)
1084: {
1085:   PC_BJacobi             *jac = (PC_BJacobi*)pc->data;
1086:   PetscErrorCode         ierr;
1087:   PetscInt               m,n_local,N,M,start,i;
1088:   const char             *prefix,*pprefix,*mprefix;
1089:   KSP                    ksp;
1090:   Vec                    x,y;
1091:   PC_BJacobi_Multiblock  *bjac = (PC_BJacobi_Multiblock*)jac->data;
1092:   PC                     subpc;
1093:   IS                     is;
1094:   MatReuse               scall = MAT_REUSE_MATRIX;

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

1099:   n_local = jac->n_local;

1101:   if (jac->use_true_local) {
1102:     PetscTruth same;
1103:     PetscTypeCompare((PetscObject)mat,((PetscObject)pmat)->type_name,&same);
1104:     if (!same) SETERRQ(PETSC_ERR_ARG_INCOMP,"Matrices not of same type");
1105:   }

1107:   if (!pc->setupcalled) {
1108:     scall                  = MAT_INITIAL_MATRIX;
1109:     pc->ops->destroy       = PCDestroy_BJacobi_Multiblock;
1110:     pc->ops->apply         = PCApply_BJacobi_Multiblock;
1111:     pc->ops->applytranspose= PCApplyTranspose_BJacobi_Multiblock;
1112:     pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiblock;

1114:     PetscNewLog(pc,PC_BJacobi_Multiblock,&bjac);
1115:     PetscMalloc(n_local*sizeof(KSP),&jac->ksp);
1116:     PetscLogObjectMemory(pc,sizeof(n_local*sizeof(KSP)));
1117:     PetscMalloc2(n_local,Vec,&bjac->x,n_local,Vec,&bjac->y);
1118:     PetscMalloc(n_local*sizeof(PetscScalar),&bjac->starts);
1119:     PetscLogObjectMemory(pc,sizeof(n_local*sizeof(PetscScalar)));
1120: 
1121:     jac->data    = (void*)bjac;
1122:     PetscMalloc(n_local*sizeof(IS),&bjac->is);
1123:     PetscLogObjectMemory(pc,sizeof(n_local*sizeof(IS)));

1125:     start = 0;
1126:     for (i=0; i<n_local; i++) {
1127:       KSPCreate(PETSC_COMM_SELF,&ksp);
1128:       PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);
1129:       PetscLogObjectParent(pc,ksp);
1130:       KSPSetType(ksp,KSPPREONLY);
1131:       KSPGetPC(ksp,&subpc);
1132:       PCGetOptionsPrefix(pc,&prefix);
1133:       KSPSetOptionsPrefix(ksp,prefix);
1134:       KSPAppendOptionsPrefix(ksp,"sub_");

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

1138:       /*
1139:       The reason we need to generate these vectors is to serve 
1140:       as the right-hand side and solution vector for the solve on the 
1141:       block. We do not need to allocate space for the vectors since
1142:       that is provided via VecPlaceArray() just before the call to 
1143:       KSPSolve() on the block.

1145:       */
1146:       VecCreateSeq(PETSC_COMM_SELF,m,&x);
1147:       VecCreateSeqWithArray(PETSC_COMM_SELF,m,PETSC_NULL,&y);
1148:       PetscLogObjectParent(pc,x);
1149:       PetscLogObjectParent(pc,y);
1150:       bjac->x[i]      = x;
1151:       bjac->y[i]      = y;
1152:       bjac->starts[i] = start;
1153:       jac->ksp[i]    = ksp;

1155:       ISCreateStride(PETSC_COMM_SELF,m,start,1,&is);
1156:       bjac->is[i] = is;
1157:       PetscLogObjectParent(pc,is);

1159:       start += m;
1160:     }
1161:   } else {
1162:     bjac = (PC_BJacobi_Multiblock*)jac->data;
1163:     /* 
1164:        Destroy the blocks from the previous iteration
1165:     */
1166:     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
1167:       MatDestroyMatrices(n_local,&bjac->pmat);
1168:       if (jac->use_true_local) {
1169:         MatDestroyMatrices(n_local,&bjac->mat);
1170:       }
1171:       scall = MAT_INITIAL_MATRIX;
1172:     }
1173:   }

1175:   MatGetSubMatrices(pmat,n_local,bjac->is,bjac->is,scall,&bjac->pmat);
1176:   if (jac->use_true_local) {
1177:     PetscObjectGetOptionsPrefix((PetscObject)mat,&mprefix);
1178:     MatGetSubMatrices(mat,n_local,bjac->is,bjac->is,scall,&bjac->mat);
1179:   }
1180:   /* Return control to the user so that the submatrices can be modified (e.g., to apply
1181:      different boundary conditions for the submatrices than for the global problem) */
1182:   PCModifySubMatrices(pc,n_local,bjac->is,bjac->is,bjac->pmat,pc->modifysubmatricesP);

1184:   PetscObjectGetOptionsPrefix((PetscObject)pmat,&pprefix);
1185:   for (i=0; i<n_local; i++) {
1186:     PetscLogObjectParent(pc,bjac->pmat[i]);
1187:     PetscObjectSetOptionsPrefix((PetscObject)bjac->pmat[i],pprefix);
1188:     if (jac->use_true_local) {
1189:       PetscLogObjectParent(pc,bjac->mat[i]);
1190:       PetscObjectSetOptionsPrefix((PetscObject)bjac->mat[i],mprefix);
1191:       KSPSetOperators(jac->ksp[i],bjac->mat[i],bjac->pmat[i],pc->flag);
1192:     } else {
1193:       KSPSetOperators(jac->ksp[i],bjac->pmat[i],bjac->pmat[i],pc->flag);
1194:     }
1195:     if(pc->setfromoptionscalled){
1196:       KSPSetFromOptions(jac->ksp[i]);
1197:     }
1198:   }
1199:   return(0);
1200: }