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: PetscOptionsName("-pc_bjacobi_truelocal","Use the true matrix, not preconditioner matrix to define matrix vector product in sub-problems","PCBJacobiSetUseTrueLocal",&flg);
220: if (flg) {
221: PCBJacobiSetUseTrueLocal(pc);
222: }
223: PetscOptionsTail();
224: return(0);
225: }
229: static PetscErrorCode PCView_BJacobi(PC pc,PetscViewer viewer)
230: {
231: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
233: PetscMPIInt rank;
234: PetscInt i;
235: PetscTruth iascii,isstring;
236: PetscViewer sviewer;
239: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
240: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
241: if (iascii) {
242: if (jac->use_true_local) {
243: PetscViewerASCIIPrintf(viewer," block Jacobi: using true local matrix, number of blocks = %D\n",jac->n);
244: }
245: PetscViewerASCIIPrintf(viewer," block Jacobi: number of blocks = %D\n",jac->n);
246: MPI_Comm_rank(((PetscObject)pc)->comm,&rank);
247: if (jac->same_local_solves) {
248: PetscViewerASCIIPrintf(viewer," Local solve is same for all blocks, in the following KSP and PC objects:\n");
249: PetscViewerGetSingleton(viewer,&sviewer);
250: if (!rank && jac->ksp) {
251: PetscViewerASCIIPushTab(viewer);
252: KSPView(jac->ksp[0],sviewer);
253: PetscViewerASCIIPopTab(viewer);
254: }
255: PetscViewerRestoreSingleton(viewer,&sviewer);
256: } else {
257: PetscInt n_global;
258: MPI_Allreduce(&jac->n_local,&n_global,1,MPIU_INT,MPI_MAX,((PetscObject)pc)->comm);
259: PetscViewerASCIIPrintf(viewer," Local solve info for each block is in the following KSP and PC objects:\n");
260: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] number of local blocks = %D, first local block number = %D\n",
261: rank,jac->n_local,jac->first_local);
262: PetscViewerASCIIPushTab(viewer);
263: for (i=0; i<n_global; i++) {
264: PetscViewerGetSingleton(viewer,&sviewer);
265: if (i < jac->n_local) {
266: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] local block number %D\n",rank,i);
267: KSPView(jac->ksp[i],sviewer);
268: PetscViewerASCIISynchronizedPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");
269: }
270: PetscViewerRestoreSingleton(viewer,&sviewer);
271: }
272: PetscViewerASCIIPopTab(viewer);
273: PetscViewerFlush(viewer);
274: }
275: } else if (isstring) {
276: PetscViewerStringSPrintf(viewer," blks=%D",jac->n);
277: PetscViewerGetSingleton(viewer,&sviewer);
278: if (jac->ksp) {KSPView(jac->ksp[0],sviewer);}
279: PetscViewerRestoreSingleton(viewer,&sviewer);
280: } else {
281: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for block Jacobi",((PetscObject)viewer)->type_name);
282: }
283: return(0);
284: }
286: /* -------------------------------------------------------------------------------------*/
291: PetscErrorCode PCBJacobiSetUseTrueLocal_BJacobi(PC pc)
292: {
293: PC_BJacobi *jac;
296: jac = (PC_BJacobi*)pc->data;
297: jac->use_true_local = PETSC_TRUE;
298: return(0);
299: }
305: PetscErrorCode PCBJacobiGetSubKSP_BJacobi(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp)
306: {
307: PC_BJacobi *jac = (PC_BJacobi*)pc->data;;
310: if (!pc->setupcalled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call KSPSetUp() or PCSetUp() first");
312: if (n_local) *n_local = jac->n_local;
313: if (first_local) *first_local = jac->first_local;
314: *ksp = jac->ksp;
315: jac->same_local_solves = PETSC_FALSE; /* Assume that local solves are now different;
316: not necessarily true though! This flag is
317: used only for PCView_BJacobi() */
318: return(0);
319: }
325: PetscErrorCode PCBJacobiSetTotalBlocks_BJacobi(PC pc,PetscInt blocks,PetscInt *lens)
326: {
327: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
332: if (pc->setupcalled > 0 && jac->n!=blocks) SETERRQ(PETSC_ERR_ORDER,"Cannot alter number of blocks after PCSetUp()/KSPSetUp() has been called");
333: jac->n = blocks;
334: if (!lens) {
335: jac->g_lens = 0;
336: } else {
337: PetscMalloc(blocks*sizeof(PetscInt),&jac->g_lens);
338: PetscLogObjectMemory(pc,blocks*sizeof(PetscInt));
339: PetscMemcpy(jac->g_lens,lens,blocks*sizeof(PetscInt));
340: }
341: return(0);
342: }
348: PetscErrorCode PCBJacobiGetTotalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[])
349: {
350: PC_BJacobi *jac = (PC_BJacobi*) pc->data;
353: *blocks = jac->n;
354: if (lens) *lens = jac->g_lens;
355: return(0);
356: }
362: PetscErrorCode PCBJacobiSetLocalBlocks_BJacobi(PC pc,PetscInt blocks,const PetscInt lens[])
363: {
364: PC_BJacobi *jac;
368: jac = (PC_BJacobi*)pc->data;
370: jac->n_local = blocks;
371: if (!lens) {
372: jac->l_lens = 0;
373: } else {
374: PetscMalloc(blocks*sizeof(PetscInt),&jac->l_lens);
375: PetscLogObjectMemory(pc,blocks*sizeof(PetscInt));
376: PetscMemcpy(jac->l_lens,lens,blocks*sizeof(PetscInt));
377: }
378: return(0);
379: }
385: PetscErrorCode PCBJacobiGetLocalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[])
386: {
387: PC_BJacobi *jac = (PC_BJacobi*) pc->data;
390: *blocks = jac->n_local;
391: if (lens) *lens = jac->l_lens;
392: return(0);
393: }
396: /* -------------------------------------------------------------------------------------*/
400: /*@
401: PCBJacobiSetUseTrueLocal - Sets a flag to indicate that the block
402: problem is associated with the linear system matrix instead of the
403: default (where it is associated with the preconditioning matrix).
404: That is, if the local system is solved iteratively then it iterates
405: on the block from the matrix using the block from the preconditioner
406: as the preconditioner for the local block.
408: Collective on PC
410: Input Parameters:
411: . pc - the preconditioner context
413: Options Database Key:
414: . -pc_bjacobi_truelocal - Activates PCBJacobiSetUseTrueLocal()
416: Notes:
417: For the common case in which the preconditioning and linear
418: system matrices are identical, this routine is unnecessary.
420: Level: intermediate
422: .keywords: block, Jacobi, set, true, local, flag
424: .seealso: PCSetOperators(), PCBJacobiSetLocalBlocks()
425: @*/
426: PetscErrorCode PCBJacobiSetUseTrueLocal(PC pc)
427: {
428: PetscErrorCode ierr,(*f)(PC);
432: PetscObjectQueryFunction((PetscObject)pc,"PCBJacobiSetUseTrueLocal_C",(void (**)(void))&f);
433: if (f) {
434: (*f)(pc);
435: }
437: return(0);
438: }
442: /*@C
443: PCBJacobiGetSubKSP - Gets the local KSP contexts for all blocks on
444: this processor.
445:
446: Note Collective
448: Input Parameter:
449: . pc - the preconditioner context
451: Output Parameters:
452: + n_local - the number of blocks on this processor, or PETSC_NULL
453: . first_local - the global number of the first block on this processor, or PETSC_NULL
454: - ksp - the array of KSP contexts
456: Notes:
457: After PCBJacobiGetSubKSP() the array of KSP contexts is not to be freed.
458:
459: Currently for some matrix implementations only 1 block per processor
460: is supported.
461:
462: You must call KSPSetUp() or PCSetUp() before calling PCBJacobiGetSubKSP().
464: Level: advanced
466: .keywords: block, Jacobi, get, sub, KSP, context
468: .seealso: PCBJacobiGetSubKSP()
469: @*/
470: PetscErrorCode PCBJacobiGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
471: {
472: PetscErrorCode ierr,(*f)(PC,PetscInt *,PetscInt *,KSP **);
476: PetscObjectQueryFunction((PetscObject)pc,"PCBJacobiGetSubKSP_C",(void (**)(void))&f);
477: if (f) {
478: (*f)(pc,n_local,first_local,ksp);
479: } else {
480: SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subsolvers for this preconditioner");
481: }
482: return(0);
483: }
487: /*@
488: PCBJacobiSetTotalBlocks - Sets the global number of blocks for the block
489: Jacobi preconditioner.
491: Collective on PC
493: Input Parameters:
494: + pc - the preconditioner context
495: . blocks - the number of blocks
496: - lens - [optional] integer array containing the size of each block
498: Options Database Key:
499: . -pc_bjacobi_blocks <blocks> - Sets the number of global blocks
501: Notes:
502: Currently only a limited number of blocking configurations are supported.
503: All processors sharing the PC must call this routine with the same data.
505: Level: intermediate
507: .keywords: set, number, Jacobi, global, total, blocks
509: .seealso: PCBJacobiSetUseTrueLocal(), PCBJacobiSetLocalBlocks()
510: @*/
511: PetscErrorCode PCBJacobiSetTotalBlocks(PC pc,PetscInt blocks,const PetscInt lens[])
512: {
513: PetscErrorCode ierr,(*f)(PC,PetscInt,const PetscInt[]);
517: if (blocks <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Must have positive blocks");
518: PetscObjectQueryFunction((PetscObject)pc,"PCBJacobiSetTotalBlocks_C",(void (**)(void))&f);
519: if (f) {
520: (*f)(pc,blocks,lens);
521: }
522: return(0);
523: }
527: /*@C
528: PCBJacobiGetTotalBlocks - Gets the global number of blocks for the block
529: Jacobi preconditioner.
531: Collective on PC
533: Input Parameter:
534: . pc - the preconditioner context
536: Output parameters:
537: + blocks - the number of blocks
538: - lens - integer array containing the size of each block
540: Level: intermediate
542: .keywords: get, number, Jacobi, global, total, blocks
544: .seealso: PCBJacobiSetUseTrueLocal(), PCBJacobiGetLocalBlocks()
545: @*/
546: PetscErrorCode PCBJacobiGetTotalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
547: {
548: PetscErrorCode ierr,(*f)(PC,PetscInt*, const PetscInt *[]);
553: PetscObjectQueryFunction((PetscObject)pc,"PCBJacobiGetTotalBlocks_C",(void (**)(void))&f);
554: if (f) {
555: (*f)(pc,blocks,lens);
556: }
557: return(0);
558: }
559:
562: /*@
563: PCBJacobiSetLocalBlocks - Sets the local number of blocks for the block
564: Jacobi preconditioner.
566: Not Collective
568: Input Parameters:
569: + pc - the preconditioner context
570: . blocks - the number of blocks
571: - lens - [optional] integer array containing size of each block
573: Note:
574: Currently only a limited number of blocking configurations are supported.
576: Level: intermediate
578: .keywords: PC, set, number, Jacobi, local, blocks
580: .seealso: PCBJacobiSetUseTrueLocal(), PCBJacobiSetTotalBlocks()
581: @*/
582: PetscErrorCode PCBJacobiSetLocalBlocks(PC pc,PetscInt blocks,const PetscInt lens[])
583: {
584: PetscErrorCode ierr,(*f)(PC,PetscInt,const PetscInt []);
588: if (blocks < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Must have nonegative blocks");
589: PetscObjectQueryFunction((PetscObject)pc,"PCBJacobiSetLocalBlocks_C",(void (**)(void))&f);
590: if (f) {
591: (*f)(pc,blocks,lens);
592: }
593: return(0);
594: }
595:
598: /*@C
599: PCBJacobiGetLocalBlocks - Gets the local number of blocks for the block
600: Jacobi preconditioner.
602: Not Collective
604: Input Parameters:
605: + pc - the preconditioner context
606: . blocks - the number of blocks
607: - lens - [optional] integer array containing size of each block
609: Note:
610: Currently only a limited number of blocking configurations are supported.
612: Level: intermediate
614: .keywords: PC, get, number, Jacobi, local, blocks
616: .seealso: PCBJacobiSetUseTrueLocal(), PCBJacobiGetTotalBlocks()
617: @*/
618: PetscErrorCode PCBJacobiGetLocalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
619: {
620: PetscErrorCode ierr,(*f)(PC,PetscInt*, const PetscInt *[]);
625: PetscObjectQueryFunction((PetscObject)pc,"PCBJacobiGetLocalBlocks_C",(void (**)(void))&f);
626: if (f) {
627: (*f)(pc,blocks,lens);
628: }
629: return(0);
630: }
632: /* -----------------------------------------------------------------------------------*/
634: /*MC
635: PCBJACOBI - Use block Jacobi preconditioning, each block is (approximately) solved with
636: its own KSP object.
638: Options Database Keys:
639: . -pc_bjacobi_truelocal - Activates PCBJacobiSetUseTrueLocal()
641: Notes: Each processor can have one or more blocks, but a block cannot be shared by more
642: than one processor. Defaults to one block per processor.
644: To set options on the solvers for each block append -sub_ to all the KSP, KSP, and PC
645: options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
646:
647: To set the options on the solvers separate for each block call PCBJacobiGetSubKSP()
648: and set the options directly on the resulting KSP object (you can access its PC
649: KSPGetPC())
651: Level: beginner
653: Concepts: block Jacobi
655: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC,
656: PCASM, PCBJacobiSetUseTrueLocal(), PCBJacobiGetSubKSP(), PCBJacobiSetTotalBlocks(),
657: PCBJacobiSetLocalBlocks(), PCSetModifySubmatrices()
658: M*/
663: PetscErrorCode PCCreate_BJacobi(PC pc)
664: {
666: PetscMPIInt rank;
667: PC_BJacobi *jac;
670: PetscNewLog(pc,PC_BJacobi,&jac);
671: MPI_Comm_rank(((PetscObject)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, 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: PetscTruth wasSetup;
895: /* set default direct solver with no Krylov method */
896: if (!pc->setupcalled) {
897: const char *prefix;
898: wasSetup = PETSC_FALSE;
899: KSPCreate(PETSC_COMM_SELF,&ksp);
900: PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);
901: PetscLogObjectParent(pc,ksp);
902: KSPSetType(ksp,KSPPREONLY);
903: PCGetOptionsPrefix(pc,&prefix);
904: KSPSetOptionsPrefix(ksp,prefix);
905: KSPAppendOptionsPrefix(ksp,"sub_");
906: /*
907: The reason we need to generate these vectors is to serve
908: as the right-hand side and solution vector for the solve on the
909: block. We do not need to allocate space for the vectors since
910: that is provided via VecPlaceArray() just before the call to
911: KSPSolve() on the block.
912: */
913: MatGetSize(pmat,&m,&m);
914: VecCreateSeqWithArray(PETSC_COMM_SELF,m,PETSC_NULL,&x);
915: VecCreateSeqWithArray(PETSC_COMM_SELF,m,PETSC_NULL,&y);
916: PetscLogObjectParent(pc,x);
917: PetscLogObjectParent(pc,y);
919: pc->ops->destroy = PCDestroy_BJacobi_Singleblock;
920: pc->ops->apply = PCApply_BJacobi_Singleblock;
921: pc->ops->applysymmetricleft = PCApplySymmetricLeft_BJacobi_Singleblock;
922: pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Singleblock;
923: pc->ops->applytranspose = PCApplyTranspose_BJacobi_Singleblock;
924: pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Singleblock;
926: PetscNewLog(pc,PC_BJacobi_Singleblock,&bjac);
927: bjac->x = x;
928: bjac->y = y;
930: PetscMalloc(sizeof(KSP),&jac->ksp);
931: jac->ksp[0] = ksp;
932: jac->data = (void*)bjac;
933: } else {
934: wasSetup = PETSC_TRUE;
935: ksp = jac->ksp[0];
936: bjac = (PC_BJacobi_Singleblock *)jac->data;
937: }
938: if (jac->use_true_local) {
939: KSPSetOperators(ksp,mat,pmat,pc->flag);
940: } else {
941: KSPSetOperators(ksp,pmat,pmat,pc->flag);
942: }
943: if (!wasSetup && pc->setfromoptionscalled) {
944: KSPSetFromOptions(ksp);
945: }
946: return(0);
947: }
949: /* ---------------------------------------------------------------------------------------------*/
953: PetscErrorCode PCDestroy_BJacobi_Multiblock(PC pc)
954: {
955: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
956: PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
957: PetscErrorCode ierr;
958: PetscInt i;
961: MatDestroyMatrices(jac->n_local,&bjac->pmat);
962: if (jac->use_true_local) {
963: MatDestroyMatrices(jac->n_local,&bjac->mat);
964: }
966: /*
967: If the on processor block had to be generated via a MatGetDiagonalBlock()
968: that creates a copy, this frees the space
969: */
970: if (jac->tp_mat) {
971: MatDestroy(jac->tp_mat);
972: }
973: if (jac->tp_pmat) {
974: MatDestroy(jac->tp_pmat);
975: }
977: for (i=0; i<jac->n_local; i++) {
978: KSPDestroy(jac->ksp[i]);
979: VecDestroy(bjac->x[i]);
980: VecDestroy(bjac->y[i]);
981: ISDestroy(bjac->is[i]);
982: }
983: PetscFree(jac->ksp);
984: PetscFree2(bjac->x,bjac->y);
985: PetscFree(bjac->starts);
986: PetscFree(bjac->is);
987: PetscFree(bjac);
988: PetscFree(jac->l_lens);
989: PetscFree(jac->g_lens);
990: PetscFree(jac);
991: return(0);
992: }
996: PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiblock(PC pc)
997: {
998: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
1000: PetscInt i,n_local = jac->n_local;
1003: for (i=0; i<n_local; i++) {
1004: KSPSetUp(jac->ksp[i]);
1005: }
1006: return(0);
1007: }
1009: /*
1010: Preconditioner for block Jacobi
1011: */
1014: PetscErrorCode PCApply_BJacobi_Multiblock(PC pc,Vec x,Vec y)
1015: {
1016: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
1017: PetscErrorCode ierr;
1018: PetscInt i,n_local = jac->n_local;
1019: PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
1020: PetscScalar *xin,*yin;
1023: VecGetArray(x,&xin);
1024: VecGetArray(y,&yin);
1025: for (i=0; i<n_local; i++) {
1026: /*
1027: To avoid copying the subvector from x into a workspace we instead
1028: make the workspace vector array point to the subpart of the array of
1029: the global vector.
1030: */
1031: VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);
1032: VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);
1034: PetscLogEventBegin(PC_SetUpOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);
1035: KSPSolve(jac->ksp[i],bjac->x[i],bjac->y[i]);
1036: PetscLogEventEnd(PC_SetUpOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);
1038: VecResetArray(bjac->x[i]);
1039: VecResetArray(bjac->y[i]);
1040: }
1041: VecRestoreArray(x,&xin);
1042: VecRestoreArray(y,&yin);
1043: return(0);
1044: }
1046: /*
1047: Preconditioner for block Jacobi
1048: */
1051: PetscErrorCode PCApplyTranspose_BJacobi_Multiblock(PC pc,Vec x,Vec y)
1052: {
1053: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
1054: PetscErrorCode ierr;
1055: PetscInt i,n_local = jac->n_local;
1056: PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
1057: PetscScalar *xin,*yin;
1060: VecGetArray(x,&xin);
1061: VecGetArray(y,&yin);
1062: for (i=0; i<n_local; i++) {
1063: /*
1064: To avoid copying the subvector from x into a workspace we instead
1065: make the workspace vector array point to the subpart of the array of
1066: the global vector.
1067: */
1068: VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);
1069: VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);
1071: PetscLogEventBegin(PC_ApplyTransposeOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);
1072: KSPSolveTranspose(jac->ksp[i],bjac->x[i],bjac->y[i]);
1073: PetscLogEventEnd(PC_ApplyTransposeOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);
1074: }
1075: VecRestoreArray(x,&xin);
1076: VecRestoreArray(y,&yin);
1077: return(0);
1078: }
1082: static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC pc,Mat mat,Mat pmat)
1083: {
1084: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
1085: PetscErrorCode ierr;
1086: PetscInt m,n_local,N,M,start,i;
1087: const char *prefix,*pprefix,*mprefix;
1088: KSP ksp;
1089: Vec x,y;
1090: PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
1091: PC subpc;
1092: IS is;
1093: MatReuse scall = MAT_REUSE_MATRIX;
1096: MatGetLocalSize(pc->pmat,&M,&N);
1098: n_local = jac->n_local;
1100: if (jac->use_true_local) {
1101: PetscTruth same;
1102: PetscTypeCompare((PetscObject)mat,((PetscObject)pmat)->type_name,&same);
1103: if (!same) SETERRQ(PETSC_ERR_ARG_INCOMP,"Matrices not of same type");
1104: }
1106: if (!pc->setupcalled) {
1107: scall = MAT_INITIAL_MATRIX;
1108: pc->ops->destroy = PCDestroy_BJacobi_Multiblock;
1109: pc->ops->apply = PCApply_BJacobi_Multiblock;
1110: pc->ops->applytranspose= PCApplyTranspose_BJacobi_Multiblock;
1111: pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiblock;
1113: PetscNewLog(pc,PC_BJacobi_Multiblock,&bjac);
1114: PetscMalloc(n_local*sizeof(KSP),&jac->ksp);
1115: PetscLogObjectMemory(pc,sizeof(n_local*sizeof(KSP)));
1116: PetscMalloc2(n_local,Vec,&bjac->x,n_local,Vec,&bjac->y);
1117: PetscMalloc(n_local*sizeof(PetscScalar),&bjac->starts);
1118: PetscLogObjectMemory(pc,sizeof(n_local*sizeof(PetscScalar)));
1119:
1120: jac->data = (void*)bjac;
1121: PetscMalloc(n_local*sizeof(IS),&bjac->is);
1122: PetscLogObjectMemory(pc,sizeof(n_local*sizeof(IS)));
1124: start = 0;
1125: for (i=0; i<n_local; i++) {
1126: KSPCreate(PETSC_COMM_SELF,&ksp);
1127: PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);
1128: PetscLogObjectParent(pc,ksp);
1129: KSPSetType(ksp,KSPPREONLY);
1130: KSPGetPC(ksp,&subpc);
1131: PCGetOptionsPrefix(pc,&prefix);
1132: KSPSetOptionsPrefix(ksp,prefix);
1133: KSPAppendOptionsPrefix(ksp,"sub_");
1135: m = jac->l_lens[i];
1137: /*
1138: The reason we need to generate these vectors is to serve
1139: as the right-hand side and solution vector for the solve on the
1140: block. We do not need to allocate space for the vectors since
1141: that is provided via VecPlaceArray() just before the call to
1142: KSPSolve() on the block.
1144: */
1145: VecCreateSeq(PETSC_COMM_SELF,m,&x);
1146: VecCreateSeqWithArray(PETSC_COMM_SELF,m,PETSC_NULL,&y);
1147: PetscLogObjectParent(pc,x);
1148: PetscLogObjectParent(pc,y);
1149: bjac->x[i] = x;
1150: bjac->y[i] = y;
1151: bjac->starts[i] = start;
1152: jac->ksp[i] = ksp;
1154: ISCreateStride(PETSC_COMM_SELF,m,start,1,&is);
1155: bjac->is[i] = is;
1156: PetscLogObjectParent(pc,is);
1158: start += m;
1159: }
1160: } else {
1161: bjac = (PC_BJacobi_Multiblock*)jac->data;
1162: /*
1163: Destroy the blocks from the previous iteration
1164: */
1165: if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
1166: MatDestroyMatrices(n_local,&bjac->pmat);
1167: if (jac->use_true_local) {
1168: MatDestroyMatrices(n_local,&bjac->mat);
1169: }
1170: scall = MAT_INITIAL_MATRIX;
1171: }
1172: }
1174: MatGetSubMatrices(pmat,n_local,bjac->is,bjac->is,scall,&bjac->pmat);
1175: if (jac->use_true_local) {
1176: PetscObjectGetOptionsPrefix((PetscObject)mat,&mprefix);
1177: MatGetSubMatrices(mat,n_local,bjac->is,bjac->is,scall,&bjac->mat);
1178: }
1179: /* Return control to the user so that the submatrices can be modified (e.g., to apply
1180: different boundary conditions for the submatrices than for the global problem) */
1181: PCModifySubMatrices(pc,n_local,bjac->is,bjac->is,bjac->pmat,pc->modifysubmatricesP);
1183: PetscObjectGetOptionsPrefix((PetscObject)pmat,&pprefix);
1184: for (i=0; i<n_local; i++) {
1185: PetscLogObjectParent(pc,bjac->pmat[i]);
1186: PetscObjectSetOptionsPrefix((PetscObject)bjac->pmat[i],pprefix);
1187: if (jac->use_true_local) {
1188: PetscLogObjectParent(pc,bjac->mat[i]);
1189: PetscObjectSetOptionsPrefix((PetscObject)bjac->mat[i],mprefix);
1190: KSPSetOperators(jac->ksp[i],bjac->mat[i],bjac->pmat[i],pc->flag);
1191: } else {
1192: KSPSetOperators(jac->ksp[i],bjac->pmat[i],bjac->pmat[i],pc->flag);
1193: }
1194: if(pc->setfromoptionscalled){
1195: KSPSetFromOptions(jac->ksp[i]);
1196: }
1197: }
1198: return(0);
1199: }