Actual source code: bjacobi.c
petsc-3.4.2 2013-07-02
2: /*
3: Defines a block Jacobi preconditioner.
4: */
5: #include <petsc-private/pcimpl.h> /*I "petscpc.h" I*/
6: #include <../src/ksp/pc/impls/bjacobi/bjacobi.h>
8: static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC,Mat,Mat);
9: static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC,Mat,Mat);
10: static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC);
14: static PetscErrorCode PCSetUp_BJacobi(PC pc)
15: {
16: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
17: Mat mat = pc->mat,pmat = pc->pmat;
18: PetscErrorCode ierr,(*f)(Mat,Mat*);
19: PetscInt N,M,start,i,sum,end;
20: PetscInt bs,i_start=-1,i_end=-1;
21: PetscMPIInt rank,size;
22: const char *pprefix,*mprefix;
25: MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);
26: MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
27: MatGetLocalSize(pc->pmat,&M,&N);
28: MatGetBlockSize(pc->pmat,&bs);
30: if (jac->n > 0 && jac->n < size) {
31: PCSetUp_BJacobi_Multiproc(pc);
32: return(0);
33: }
35: /* --------------------------------------------------------------------------
36: Determines the number of blocks assigned to each processor
37: -----------------------------------------------------------------------------*/
39: /* local block count given */
40: if (jac->n_local > 0 && jac->n < 0) {
41: MPI_Allreduce(&jac->n_local,&jac->n,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));
42: if (jac->l_lens) { /* check that user set these correctly */
43: sum = 0;
44: for (i=0; i<jac->n_local; i++) {
45: if (jac->l_lens[i]/bs*bs !=jac->l_lens[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat blocksize doesn't match block Jacobi layout");
46: sum += jac->l_lens[i];
47: }
48: if (sum != M) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local lens set incorrectly");
49: } else {
50: PetscMalloc(jac->n_local*sizeof(PetscInt),&jac->l_lens);
51: for (i=0; i<jac->n_local; i++) jac->l_lens[i] = bs*((M/bs)/jac->n_local + (((M/bs) % jac->n_local) > i));
52: }
53: } else if (jac->n > 0 && jac->n_local < 0) { /* global block count given */
54: /* global blocks given: determine which ones are local */
55: if (jac->g_lens) {
56: /* check if the g_lens is has valid entries */
57: for (i=0; i<jac->n; i++) {
58: if (!jac->g_lens[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Zero block not allowed");
59: if (jac->g_lens[i]/bs*bs != jac->g_lens[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat blocksize doesn't match block Jacobi layout");
60: }
61: if (size == 1) {
62: jac->n_local = jac->n;
63: PetscMalloc(jac->n_local*sizeof(PetscInt),&jac->l_lens);
64: PetscMemcpy(jac->l_lens,jac->g_lens,jac->n_local*sizeof(PetscInt));
65: /* check that user set these correctly */
66: sum = 0;
67: for (i=0; i<jac->n_local; i++) sum += jac->l_lens[i];
68: if (sum != M) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Global lens set incorrectly");
69: } else {
70: MatGetOwnershipRange(pc->pmat,&start,&end);
71: /* loop over blocks determing first one owned by me */
72: sum = 0;
73: for (i=0; i<jac->n+1; i++) {
74: if (sum == start) { i_start = i; goto start_1;}
75: if (i < jac->n) sum += jac->g_lens[i];
76: }
77: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Block sizes used in PCBJacobiSetTotalBlocks()\nare not compatible with parallel matrix layout");
78: start_1:
79: for (i=i_start; i<jac->n+1; i++) {
80: if (sum == end) { i_end = i; goto end_1; }
81: if (i < jac->n) sum += jac->g_lens[i];
82: }
83: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Block sizes used in PCBJacobiSetTotalBlocks()\nare not compatible with parallel matrix layout");
84: end_1:
85: jac->n_local = i_end - i_start;
86: PetscMalloc(jac->n_local*sizeof(PetscInt),&jac->l_lens);
87: PetscMemcpy(jac->l_lens,jac->g_lens+i_start,jac->n_local*sizeof(PetscInt));
88: }
89: } else { /* no global blocks given, determine then using default layout */
90: jac->n_local = jac->n/size + ((jac->n % size) > rank);
91: PetscMalloc(jac->n_local*sizeof(PetscInt),&jac->l_lens);
92: for (i=0; i<jac->n_local; i++) {
93: jac->l_lens[i] = ((M/bs)/jac->n_local + (((M/bs) % jac->n_local) > i))*bs;
94: if (!jac->l_lens[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Too many blocks given");
95: }
96: }
97: } else if (jac->n < 0 && jac->n_local < 0) { /* no blocks given */
98: jac->n = size;
99: jac->n_local = 1;
100: PetscMalloc(sizeof(PetscInt),&jac->l_lens);
101: jac->l_lens[0] = M;
102: } else { /* jac->n > 0 && jac->n_local > 0 */
103: if (!jac->l_lens) {
104: PetscMalloc(jac->n_local*sizeof(PetscInt),&jac->l_lens);
105: for (i=0; i<jac->n_local; i++) jac->l_lens[i] = bs*((M/bs)/jac->n_local + (((M/bs) % jac->n_local) > i));
106: }
107: }
108: if (jac->n_local < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of blocks is less than number of processors");
110: /* -------------------------
111: Determines mat and pmat
112: ---------------------------*/
113: PetscObjectQueryFunction((PetscObject)pc->mat,"MatGetDiagonalBlock_C",&f);
114: if (!f && size == 1) {
115: mat = pc->mat;
116: pmat = pc->pmat;
117: } else {
118: if (pc->useAmat) {
119: /* use block from Amat matrix, not Pmat for local MatMult() */
120: MatGetDiagonalBlock(pc->mat,&mat);
121: /* make submatrix have same prefix as entire matrix */
122: PetscObjectGetOptionsPrefix((PetscObject)pc->mat,&mprefix);
123: PetscObjectSetOptionsPrefix((PetscObject)mat,mprefix);
124: }
125: if (pc->pmat != pc->mat || !pc->useAmat) {
126: MatGetDiagonalBlock(pc->pmat,&pmat);
127: /* make submatrix have same prefix as entire matrix */
128: PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);
129: PetscObjectSetOptionsPrefix((PetscObject)pmat,pprefix);
130: } else pmat = mat;
131: }
133: /* ------
134: Setup code depends on the number of blocks
135: */
136: if (jac->n_local == 1) {
137: PCSetUp_BJacobi_Singleblock(pc,mat,pmat);
138: } else {
139: PCSetUp_BJacobi_Multiblock(pc,mat,pmat);
140: }
141: return(0);
142: }
144: /* Default destroy, if it has never been setup */
147: static PetscErrorCode PCDestroy_BJacobi(PC pc)
148: {
149: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
153: PetscFree(jac->g_lens);
154: PetscFree(jac->l_lens);
155: PetscFree(pc->data);
156: return(0);
157: }
162: static PetscErrorCode PCSetFromOptions_BJacobi(PC pc)
163: {
164: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
166: PetscInt blocks;
167: PetscBool flg;
170: PetscOptionsHead("Block Jacobi options");
171: PetscOptionsInt("-pc_bjacobi_blocks","Total number of blocks","PCBJacobiSetTotalBlocks",jac->n,&blocks,&flg);
172: if (flg) {
173: PCBJacobiSetTotalBlocks(pc,blocks,NULL);
174: }
175: PetscOptionsTail();
176: return(0);
177: }
179: #include <petscdraw.h>
182: static PetscErrorCode PCView_BJacobi(PC pc,PetscViewer viewer)
183: {
184: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
185: PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
186: PetscErrorCode ierr;
187: PetscMPIInt rank;
188: PetscInt i;
189: PetscBool iascii,isstring,isdraw;
190: PetscViewer sviewer;
193: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
194: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
195: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
196: if (iascii) {
197: if (pc->useAmat) {
198: PetscViewerASCIIPrintf(viewer," block Jacobi: using Amat local matrix, number of blocks = %D\n",jac->n);
199: }
200: PetscViewerASCIIPrintf(viewer," block Jacobi: number of blocks = %D\n",jac->n);
201: MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);
202: if (jac->same_local_solves) {
203: PetscViewerASCIIPrintf(viewer," Local solve is same for all blocks, in the following KSP and PC objects:\n");
204: if (jac->ksp && !jac->psubcomm) {
205: PetscViewerGetSingleton(viewer,&sviewer);
206: if (!rank) {
207: PetscViewerASCIIPushTab(viewer);
208: KSPView(jac->ksp[0],sviewer);
209: PetscViewerASCIIPopTab(viewer);
210: }
211: PetscViewerRestoreSingleton(viewer,&sviewer);
212: } else if (jac->psubcomm && !jac->psubcomm->color) {
213: PetscViewerASCIIGetStdout(mpjac->psubcomm->comm,&sviewer);
214: PetscViewerASCIIPushTab(viewer);
215: KSPView(*(jac->ksp),sviewer);
216: PetscViewerASCIIPopTab(viewer);
217: }
218: } else {
219: PetscInt n_global;
220: MPI_Allreduce(&jac->n_local,&n_global,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)pc));
221: PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);
222: PetscViewerASCIIPrintf(viewer," Local solve info for each block is in the following KSP and PC objects:\n");
223: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] number of local blocks = %D, first local block number = %D\n",
224: rank,jac->n_local,jac->first_local);
225: PetscViewerASCIIPushTab(viewer);
226: for (i=0; i<n_global; i++) {
227: PetscViewerGetSingleton(viewer,&sviewer);
228: if (i < jac->n_local) {
229: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] local block number %D\n",rank,i);
230: KSPView(jac->ksp[i],sviewer);
231: PetscViewerASCIISynchronizedPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");
232: }
233: PetscViewerRestoreSingleton(viewer,&sviewer);
234: }
235: PetscViewerASCIIPopTab(viewer);
236: PetscViewerFlush(viewer);
237: PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);
238: }
239: } else if (isstring) {
240: PetscViewerStringSPrintf(viewer," blks=%D",jac->n);
241: PetscViewerGetSingleton(viewer,&sviewer);
242: if (jac->ksp) {KSPView(jac->ksp[0],sviewer);}
243: PetscViewerRestoreSingleton(viewer,&sviewer);
244: } else if (isdraw) {
245: PetscDraw draw;
246: char str[25];
247: PetscReal x,y,bottom,h;
249: PetscViewerDrawGetDraw(viewer,0,&draw);
250: PetscDrawGetCurrentPoint(draw,&x,&y);
251: PetscSNPrintf(str,25,"Number blocks %D",jac->n);
252: PetscDrawBoxedString(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);
253: bottom = y - h;
254: PetscDrawPushCurrentPoint(draw,x,bottom);
255: /* warning the communicator on viewer is different then on ksp in parallel */
256: if (jac->ksp) {KSPView(jac->ksp[0],viewer);}
257: PetscDrawPopCurrentPoint(draw);
258: }
259: return(0);
260: }
262: /* -------------------------------------------------------------------------------------*/
266: static PetscErrorCode PCBJacobiGetSubKSP_BJacobi(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp)
267: {
268: PC_BJacobi *jac = (PC_BJacobi*)pc->data;;
271: if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call KSPSetUp() or PCSetUp() first");
273: if (n_local) *n_local = jac->n_local;
274: if (first_local) *first_local = jac->first_local;
275: *ksp = jac->ksp;
276: jac->same_local_solves = PETSC_FALSE; /* Assume that local solves are now different;
277: not necessarily true though! This flag is
278: used only for PCView_BJacobi() */
279: return(0);
280: }
284: static PetscErrorCode PCBJacobiSetTotalBlocks_BJacobi(PC pc,PetscInt blocks,PetscInt *lens)
285: {
286: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
290: if (pc->setupcalled > 0 && jac->n!=blocks) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Cannot alter number of blocks after PCSetUp()/KSPSetUp() has been called");
291: jac->n = blocks;
292: if (!lens) jac->g_lens = 0;
293: else {
294: PetscMalloc(blocks*sizeof(PetscInt),&jac->g_lens);
295: PetscLogObjectMemory(pc,blocks*sizeof(PetscInt));
296: PetscMemcpy(jac->g_lens,lens,blocks*sizeof(PetscInt));
297: }
298: return(0);
299: }
303: static PetscErrorCode PCBJacobiGetTotalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[])
304: {
305: PC_BJacobi *jac = (PC_BJacobi*) pc->data;
308: *blocks = jac->n;
309: if (lens) *lens = jac->g_lens;
310: return(0);
311: }
315: static PetscErrorCode PCBJacobiSetLocalBlocks_BJacobi(PC pc,PetscInt blocks,const PetscInt lens[])
316: {
317: PC_BJacobi *jac;
321: jac = (PC_BJacobi*)pc->data;
323: jac->n_local = blocks;
324: if (!lens) jac->l_lens = 0;
325: else {
326: PetscMalloc(blocks*sizeof(PetscInt),&jac->l_lens);
327: PetscLogObjectMemory(pc,blocks*sizeof(PetscInt));
328: PetscMemcpy(jac->l_lens,lens,blocks*sizeof(PetscInt));
329: }
330: return(0);
331: }
335: static PetscErrorCode PCBJacobiGetLocalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[])
336: {
337: PC_BJacobi *jac = (PC_BJacobi*) pc->data;
340: *blocks = jac->n_local;
341: if (lens) *lens = jac->l_lens;
342: return(0);
343: }
345: /* -------------------------------------------------------------------------------------*/
349: /*@C
350: PCBJacobiGetSubKSP - Gets the local KSP contexts for all blocks on
351: this processor.
353: Note Collective
355: Input Parameter:
356: . pc - the preconditioner context
358: Output Parameters:
359: + n_local - the number of blocks on this processor, or NULL
360: . first_local - the global number of the first block on this processor, or NULL
361: - ksp - the array of KSP contexts
363: Notes:
364: After PCBJacobiGetSubKSP() the array of KSP contexts is not to be freed.
366: Currently for some matrix implementations only 1 block per processor
367: is supported.
369: You must call KSPSetUp() or PCSetUp() before calling PCBJacobiGetSubKSP().
371: Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs.
372: You can call PCBJacobiGetSubKSP(pc,nlocal,firstlocal,NULL_OBJECT,ierr) to determine how large the
373: KSP array must be.
375: Level: advanced
377: .keywords: block, Jacobi, get, sub, KSP, context
379: .seealso: PCBJacobiGetSubKSP()
380: @*/
381: PetscErrorCode PCBJacobiGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
382: {
387: PetscUseMethod(pc,"PCBJacobiGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));
388: return(0);
389: }
393: /*@
394: PCBJacobiSetTotalBlocks - Sets the global number of blocks for the block
395: Jacobi preconditioner.
397: Collective on PC
399: Input Parameters:
400: + pc - the preconditioner context
401: . blocks - the number of blocks
402: - lens - [optional] integer array containing the size of each block
404: Options Database Key:
405: . -pc_bjacobi_blocks <blocks> - Sets the number of global blocks
407: Notes:
408: Currently only a limited number of blocking configurations are supported.
409: All processors sharing the PC must call this routine with the same data.
411: Level: intermediate
413: .keywords: set, number, Jacobi, global, total, blocks
415: .seealso: PCSetUseAmat(), PCBJacobiSetLocalBlocks()
416: @*/
417: PetscErrorCode PCBJacobiSetTotalBlocks(PC pc,PetscInt blocks,const PetscInt lens[])
418: {
423: if (blocks <= 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Must have positive blocks");
424: PetscTryMethod(pc,"PCBJacobiSetTotalBlocks_C",(PC,PetscInt,const PetscInt[]),(pc,blocks,lens));
425: return(0);
426: }
430: /*@C
431: PCBJacobiGetTotalBlocks - Gets the global number of blocks for the block
432: Jacobi preconditioner.
434: Not Collective
436: Input Parameter:
437: . pc - the preconditioner context
439: Output parameters:
440: + blocks - the number of blocks
441: - lens - integer array containing the size of each block
443: Level: intermediate
445: .keywords: get, number, Jacobi, global, total, blocks
447: .seealso: PCSetUseAmat(), PCBJacobiGetLocalBlocks()
448: @*/
449: PetscErrorCode PCBJacobiGetTotalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
450: {
456: PetscUseMethod(pc,"PCBJacobiGetTotalBlocks_C",(PC,PetscInt*, const PetscInt *[]),(pc,blocks,lens));
457: return(0);
458: }
462: /*@
463: PCBJacobiSetLocalBlocks - Sets the local number of blocks for the block
464: Jacobi preconditioner.
466: Not Collective
468: Input Parameters:
469: + pc - the preconditioner context
470: . blocks - the number of blocks
471: - lens - [optional] integer array containing size of each block
473: Note:
474: Currently only a limited number of blocking configurations are supported.
476: Level: intermediate
478: .keywords: PC, set, number, Jacobi, local, blocks
480: .seealso: PCSetUseAmat(), PCBJacobiSetTotalBlocks()
481: @*/
482: PetscErrorCode PCBJacobiSetLocalBlocks(PC pc,PetscInt blocks,const PetscInt lens[])
483: {
488: if (blocks < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have nonegative blocks");
489: PetscTryMethod(pc,"PCBJacobiSetLocalBlocks_C",(PC,PetscInt,const PetscInt []),(pc,blocks,lens));
490: return(0);
491: }
495: /*@C
496: PCBJacobiGetLocalBlocks - Gets the local number of blocks for the block
497: Jacobi preconditioner.
499: Not Collective
501: Input Parameters:
502: + pc - the preconditioner context
503: . blocks - the number of blocks
504: - lens - [optional] integer array containing size of each block
506: Note:
507: Currently only a limited number of blocking configurations are supported.
509: Level: intermediate
511: .keywords: PC, get, number, Jacobi, local, blocks
513: .seealso: PCSetUseAmat(), PCBJacobiGetTotalBlocks()
514: @*/
515: PetscErrorCode PCBJacobiGetLocalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
516: {
522: PetscUseMethod(pc,"PCBJacobiGetLocalBlocks_C",(PC,PetscInt*, const PetscInt *[]),(pc,blocks,lens));
523: return(0);
524: }
526: /* -----------------------------------------------------------------------------------*/
528: /*MC
529: PCBJACOBI - Use block Jacobi preconditioning, each block is (approximately) solved with
530: its own KSP object.
532: Options Database Keys:
533: . -pc_use_amat - use Amat to apply block of operator in inner Krylov method
535: Notes: Each processor can have one or more blocks, but a block cannot be shared by more
536: than one processor. Defaults to one block per processor.
538: To set options on the solvers for each block append -sub_ to all the KSP, KSP, and PC
539: options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
541: To set the options on the solvers separate for each block call PCBJacobiGetSubKSP()
542: and set the options directly on the resulting KSP object (you can access its PC
543: KSPGetPC())
545: Level: beginner
547: Concepts: block Jacobi
549: Developer Notes: This preconditioner does not currently work with CUDA/CUSP for a couple of reasons.
550: (1) It creates seq vectors as work vectors that should be cusp
551: (2) The use of VecPlaceArray() is not handled properly by CUSP (that is it will not know where
552: the ownership of the vector is so may use wrong values) even if it did know the ownership
553: it may induce extra copy ups and downs. Satish suggests a VecTransplantArray() to handle two
554: vectors sharing the same pointer and handling the CUSP side as well instead of VecGetArray()/VecPlaceArray().
557: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC,
558: PCASM, PCSetUseAmat(), PCGetUseAmat(), PCBJacobiGetSubKSP(), PCBJacobiSetTotalBlocks(),
559: PCBJacobiSetLocalBlocks(), PCSetModifySubmatrices()
560: M*/
564: PETSC_EXTERN PetscErrorCode PCCreate_BJacobi(PC pc)
565: {
567: PetscMPIInt rank;
568: PC_BJacobi *jac;
571: PetscNewLog(pc,PC_BJacobi,&jac);
572: MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);
574: pc->ops->apply = 0;
575: pc->ops->applytranspose = 0;
576: pc->ops->setup = PCSetUp_BJacobi;
577: pc->ops->destroy = PCDestroy_BJacobi;
578: pc->ops->setfromoptions = PCSetFromOptions_BJacobi;
579: pc->ops->view = PCView_BJacobi;
580: pc->ops->applyrichardson = 0;
582: pc->data = (void*)jac;
583: jac->n = -1;
584: jac->n_local = -1;
585: jac->first_local = rank;
586: jac->ksp = 0;
587: jac->same_local_solves = PETSC_TRUE;
588: jac->g_lens = 0;
589: jac->l_lens = 0;
590: jac->psubcomm = 0;
592: PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiGetSubKSP_C",PCBJacobiGetSubKSP_BJacobi);
593: PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiSetTotalBlocks_C",PCBJacobiSetTotalBlocks_BJacobi);
594: PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiGetTotalBlocks_C",PCBJacobiGetTotalBlocks_BJacobi);
595: PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiSetLocalBlocks_C",PCBJacobiSetLocalBlocks_BJacobi);
596: PetscObjectComposeFunction((PetscObject)pc,"PCBJacobiGetLocalBlocks_C",PCBJacobiGetLocalBlocks_BJacobi);
597: return(0);
598: }
600: /* --------------------------------------------------------------------------------------------*/
601: /*
602: These are for a single block per processor; works for AIJ, BAIJ; Seq and MPI
603: */
606: PetscErrorCode PCReset_BJacobi_Singleblock(PC pc)
607: {
608: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
609: PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
610: PetscErrorCode ierr;
613: KSPReset(jac->ksp[0]);
614: VecDestroy(&bjac->x);
615: VecDestroy(&bjac->y);
616: return(0);
617: }
621: PetscErrorCode PCDestroy_BJacobi_Singleblock(PC pc)
622: {
623: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
624: PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
625: PetscErrorCode ierr;
628: PCReset_BJacobi_Singleblock(pc);
629: KSPDestroy(&jac->ksp[0]);
630: PetscFree(jac->ksp);
631: PetscFree(jac->l_lens);
632: PetscFree(jac->g_lens);
633: PetscFree(bjac);
634: PetscFree(pc->data);
635: return(0);
636: }
640: PetscErrorCode PCSetUpOnBlocks_BJacobi_Singleblock(PC pc)
641: {
643: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
646: KSPSetUp(jac->ksp[0]);
647: return(0);
648: }
652: PetscErrorCode PCApply_BJacobi_Singleblock(PC pc,Vec x,Vec y)
653: {
654: PetscErrorCode ierr;
655: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
656: PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
657: PetscScalar *x_array,*y_array;
660: /*
661: The VecPlaceArray() is to avoid having to copy the
662: y vector into the bjac->x vector. The reason for
663: the bjac->x vector is that we need a sequential vector
664: for the sequential solve.
665: */
666: VecGetArray(x,&x_array);
667: VecGetArray(y,&y_array);
668: VecPlaceArray(bjac->x,x_array);
669: VecPlaceArray(bjac->y,y_array);
670: KSPSolve(jac->ksp[0],bjac->x,bjac->y);
671: VecResetArray(bjac->x);
672: VecResetArray(bjac->y);
673: VecRestoreArray(x,&x_array);
674: VecRestoreArray(y,&y_array);
675: return(0);
676: }
680: PetscErrorCode PCApplySymmetricLeft_BJacobi_Singleblock(PC pc,Vec x,Vec y)
681: {
682: PetscErrorCode ierr;
683: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
684: PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
685: PetscScalar *x_array,*y_array;
686: PC subpc;
689: /*
690: The VecPlaceArray() is to avoid having to copy the
691: y vector into the bjac->x vector. The reason for
692: the bjac->x vector is that we need a sequential vector
693: for the sequential solve.
694: */
695: VecGetArray(x,&x_array);
696: VecGetArray(y,&y_array);
697: VecPlaceArray(bjac->x,x_array);
698: VecPlaceArray(bjac->y,y_array);
699: /* apply the symmetric left portion of the inner PC operator */
700: /* note this by-passes the inner KSP and its options completely */
701: KSPGetPC(jac->ksp[0],&subpc);
702: PCApplySymmetricLeft(subpc,bjac->x,bjac->y);
703: VecResetArray(bjac->x);
704: VecResetArray(bjac->y);
705: VecRestoreArray(x,&x_array);
706: VecRestoreArray(y,&y_array);
707: return(0);
708: }
712: PetscErrorCode PCApplySymmetricRight_BJacobi_Singleblock(PC pc,Vec x,Vec y)
713: {
714: PetscErrorCode ierr;
715: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
716: PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
717: PetscScalar *x_array,*y_array;
718: PC subpc;
721: /*
722: The VecPlaceArray() is to avoid having to copy the
723: y vector into the bjac->x vector. The reason for
724: the bjac->x vector is that we need a sequential vector
725: for the sequential solve.
726: */
727: VecGetArray(x,&x_array);
728: VecGetArray(y,&y_array);
729: VecPlaceArray(bjac->x,x_array);
730: VecPlaceArray(bjac->y,y_array);
732: /* apply the symmetric right portion of the inner PC operator */
733: /* note this by-passes the inner KSP and its options completely */
735: KSPGetPC(jac->ksp[0],&subpc);
736: PCApplySymmetricRight(subpc,bjac->x,bjac->y);
738: VecRestoreArray(x,&x_array);
739: VecRestoreArray(y,&y_array);
740: return(0);
741: }
745: PetscErrorCode PCApplyTranspose_BJacobi_Singleblock(PC pc,Vec x,Vec y)
746: {
747: PetscErrorCode ierr;
748: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
749: PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data;
750: PetscScalar *x_array,*y_array;
753: /*
754: The VecPlaceArray() is to avoid having to copy the
755: y vector into the bjac->x vector. The reason for
756: the bjac->x vector is that we need a sequential vector
757: for the sequential solve.
758: */
759: VecGetArray(x,&x_array);
760: VecGetArray(y,&y_array);
761: VecPlaceArray(bjac->x,x_array);
762: VecPlaceArray(bjac->y,y_array);
763: KSPSolveTranspose(jac->ksp[0],bjac->x,bjac->y);
764: VecResetArray(bjac->x);
765: VecResetArray(bjac->y);
766: VecRestoreArray(x,&x_array);
767: VecRestoreArray(y,&y_array);
768: return(0);
769: }
773: static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC pc,Mat mat,Mat pmat)
774: {
775: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
776: PetscErrorCode ierr;
777: PetscInt m;
778: KSP ksp;
779: PC_BJacobi_Singleblock *bjac;
780: PetscBool wasSetup = PETSC_TRUE;
783: if (!pc->setupcalled) {
784: const char *prefix;
786: if (!jac->ksp) {
787: wasSetup = PETSC_FALSE;
789: KSPCreate(PETSC_COMM_SELF,&ksp);
790: PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);
791: PetscLogObjectParent(pc,ksp);
792: KSPSetType(ksp,KSPPREONLY);
793: PCGetOptionsPrefix(pc,&prefix);
794: KSPSetOptionsPrefix(ksp,prefix);
795: KSPAppendOptionsPrefix(ksp,"sub_");
797: pc->ops->reset = PCReset_BJacobi_Singleblock;
798: pc->ops->destroy = PCDestroy_BJacobi_Singleblock;
799: pc->ops->apply = PCApply_BJacobi_Singleblock;
800: pc->ops->applysymmetricleft = PCApplySymmetricLeft_BJacobi_Singleblock;
801: pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Singleblock;
802: pc->ops->applytranspose = PCApplyTranspose_BJacobi_Singleblock;
803: pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Singleblock;
805: PetscMalloc(sizeof(KSP),&jac->ksp);
806: jac->ksp[0] = ksp;
808: PetscNewLog(pc,PC_BJacobi_Singleblock,&bjac);
809: jac->data = (void*)bjac;
810: } else {
811: ksp = jac->ksp[0];
812: bjac = (PC_BJacobi_Singleblock*)jac->data;
813: }
815: /*
816: The reason we need to generate these vectors is to serve
817: as the right-hand side and solution vector for the solve on the
818: block. We do not need to allocate space for the vectors since
819: that is provided via VecPlaceArray() just before the call to
820: KSPSolve() on the block.
821: */
822: MatGetSize(pmat,&m,&m);
823: VecCreateSeqWithArray(PETSC_COMM_SELF,1,m,NULL,&bjac->x);
824: VecCreateSeqWithArray(PETSC_COMM_SELF,1,m,NULL,&bjac->y);
825: PetscLogObjectParent(pc,bjac->x);
826: PetscLogObjectParent(pc,bjac->y);
827: } else {
828: ksp = jac->ksp[0];
829: bjac = (PC_BJacobi_Singleblock*)jac->data;
830: }
831: if (pc->useAmat) {
832: KSPSetOperators(ksp,mat,pmat,pc->flag);
833: } else {
834: KSPSetOperators(ksp,pmat,pmat,pc->flag);
835: }
836: if (!wasSetup && pc->setfromoptionscalled) {
837: KSPSetFromOptions(ksp);
838: }
839: return(0);
840: }
842: /* ---------------------------------------------------------------------------------------------*/
845: PetscErrorCode PCReset_BJacobi_Multiblock(PC pc)
846: {
847: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
848: PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
849: PetscErrorCode ierr;
850: PetscInt i;
853: if (bjac && bjac->pmat) {
854: MatDestroyMatrices(jac->n_local,&bjac->pmat);
855: if (pc->useAmat) {
856: MatDestroyMatrices(jac->n_local,&bjac->mat);
857: }
858: }
860: for (i=0; i<jac->n_local; i++) {
861: KSPReset(jac->ksp[i]);
862: if (bjac && bjac->x) {
863: VecDestroy(&bjac->x[i]);
864: VecDestroy(&bjac->y[i]);
865: ISDestroy(&bjac->is[i]);
866: }
867: }
868: PetscFree(jac->l_lens);
869: PetscFree(jac->g_lens);
870: return(0);
871: }
875: PetscErrorCode PCDestroy_BJacobi_Multiblock(PC pc)
876: {
877: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
878: PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
879: PetscErrorCode ierr;
880: PetscInt i;
883: PCReset_BJacobi_Multiblock(pc);
884: if (bjac) {
885: PetscFree2(bjac->x,bjac->y);
886: PetscFree(bjac->starts);
887: PetscFree(bjac->is);
888: }
889: PetscFree(jac->data);
890: for (i=0; i<jac->n_local; i++) {
891: KSPDestroy(&jac->ksp[i]);
892: }
893: PetscFree(jac->ksp);
894: PetscFree(pc->data);
895: return(0);
896: }
900: PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiblock(PC pc)
901: {
902: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
904: PetscInt i,n_local = jac->n_local;
907: for (i=0; i<n_local; i++) {
908: KSPSetUp(jac->ksp[i]);
909: }
910: return(0);
911: }
913: /*
914: Preconditioner for block Jacobi
915: */
918: PetscErrorCode PCApply_BJacobi_Multiblock(PC pc,Vec x,Vec y)
919: {
920: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
921: PetscErrorCode ierr;
922: PetscInt i,n_local = jac->n_local;
923: PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
924: PetscScalar *xin,*yin;
927: VecGetArray(x,&xin);
928: VecGetArray(y,&yin);
929: for (i=0; i<n_local; i++) {
930: /*
931: To avoid copying the subvector from x into a workspace we instead
932: make the workspace vector array point to the subpart of the array of
933: the global vector.
934: */
935: VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);
936: VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);
938: PetscLogEventBegin(PC_ApplyOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);
939: KSPSolve(jac->ksp[i],bjac->x[i],bjac->y[i]);
940: PetscLogEventEnd(PC_ApplyOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);
942: VecResetArray(bjac->x[i]);
943: VecResetArray(bjac->y[i]);
944: }
945: VecRestoreArray(x,&xin);
946: VecRestoreArray(y,&yin);
947: return(0);
948: }
950: /*
951: Preconditioner for block Jacobi
952: */
955: PetscErrorCode PCApplyTranspose_BJacobi_Multiblock(PC pc,Vec x,Vec y)
956: {
957: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
958: PetscErrorCode ierr;
959: PetscInt i,n_local = jac->n_local;
960: PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
961: PetscScalar *xin,*yin;
964: VecGetArray(x,&xin);
965: VecGetArray(y,&yin);
966: for (i=0; i<n_local; i++) {
967: /*
968: To avoid copying the subvector from x into a workspace we instead
969: make the workspace vector array point to the subpart of the array of
970: the global vector.
971: */
972: VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);
973: VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);
975: PetscLogEventBegin(PC_ApplyTransposeOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);
976: KSPSolveTranspose(jac->ksp[i],bjac->x[i],bjac->y[i]);
977: PetscLogEventEnd(PC_ApplyTransposeOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);
979: VecResetArray(bjac->x[i]);
980: VecResetArray(bjac->y[i]);
981: }
982: VecRestoreArray(x,&xin);
983: VecRestoreArray(y,&yin);
984: return(0);
985: }
989: static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC pc,Mat mat,Mat pmat)
990: {
991: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
992: PetscErrorCode ierr;
993: PetscInt m,n_local,N,M,start,i;
994: const char *prefix,*pprefix,*mprefix;
995: KSP ksp;
996: Vec x,y;
997: PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data;
998: PC subpc;
999: IS is;
1000: MatReuse scall;
1003: MatGetLocalSize(pc->pmat,&M,&N);
1005: n_local = jac->n_local;
1007: if (pc->useAmat) {
1008: PetscBool same;
1009: PetscObjectTypeCompare((PetscObject)mat,((PetscObject)pmat)->type_name,&same);
1010: if (!same) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Matrices not of same type");
1011: }
1013: if (!pc->setupcalled) {
1014: scall = MAT_INITIAL_MATRIX;
1016: if (!jac->ksp) {
1017: pc->ops->reset = PCReset_BJacobi_Multiblock;
1018: pc->ops->destroy = PCDestroy_BJacobi_Multiblock;
1019: pc->ops->apply = PCApply_BJacobi_Multiblock;
1020: pc->ops->applytranspose= PCApplyTranspose_BJacobi_Multiblock;
1021: pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiblock;
1023: PetscNewLog(pc,PC_BJacobi_Multiblock,&bjac);
1024: PetscMalloc(n_local*sizeof(KSP),&jac->ksp);
1025: PetscLogObjectMemory(pc,sizeof(n_local*sizeof(KSP)));
1026: PetscMalloc2(n_local,Vec,&bjac->x,n_local,Vec,&bjac->y);
1027: PetscMalloc(n_local*sizeof(PetscScalar),&bjac->starts);
1028: PetscLogObjectMemory(pc,sizeof(n_local*sizeof(PetscScalar)));
1030: jac->data = (void*)bjac;
1031: PetscMalloc(n_local*sizeof(IS),&bjac->is);
1032: PetscLogObjectMemory(pc,sizeof(n_local*sizeof(IS)));
1034: for (i=0; i<n_local; i++) {
1035: KSPCreate(PETSC_COMM_SELF,&ksp);
1036: PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);
1037: PetscLogObjectParent(pc,ksp);
1038: KSPSetType(ksp,KSPPREONLY);
1039: KSPGetPC(ksp,&subpc);
1040: PCGetOptionsPrefix(pc,&prefix);
1041: KSPSetOptionsPrefix(ksp,prefix);
1042: KSPAppendOptionsPrefix(ksp,"sub_");
1044: jac->ksp[i] = ksp;
1045: }
1046: } else {
1047: bjac = (PC_BJacobi_Multiblock*)jac->data;
1048: }
1050: start = 0;
1051: for (i=0; i<n_local; i++) {
1052: m = jac->l_lens[i];
1053: /*
1054: The reason we need to generate these vectors is to serve
1055: as the right-hand side and solution vector for the solve on the
1056: block. We do not need to allocate space for the vectors since
1057: that is provided via VecPlaceArray() just before the call to
1058: KSPSolve() on the block.
1060: */
1061: VecCreateSeq(PETSC_COMM_SELF,m,&x);
1062: VecCreateSeqWithArray(PETSC_COMM_SELF,1,m,NULL,&y);
1063: PetscLogObjectParent(pc,x);
1064: PetscLogObjectParent(pc,y);
1066: bjac->x[i] = x;
1067: bjac->y[i] = y;
1068: bjac->starts[i] = start;
1070: ISCreateStride(PETSC_COMM_SELF,m,start,1,&is);
1071: bjac->is[i] = is;
1072: PetscLogObjectParent(pc,is);
1074: start += m;
1075: }
1076: } else {
1077: bjac = (PC_BJacobi_Multiblock*)jac->data;
1078: /*
1079: Destroy the blocks from the previous iteration
1080: */
1081: if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
1082: MatDestroyMatrices(n_local,&bjac->pmat);
1083: if (pc->useAmat) {
1084: MatDestroyMatrices(n_local,&bjac->mat);
1085: }
1086: scall = MAT_INITIAL_MATRIX;
1087: } else scall = MAT_REUSE_MATRIX;
1088: }
1090: MatGetSubMatrices(pmat,n_local,bjac->is,bjac->is,scall,&bjac->pmat);
1091: if (pc->useAmat) {
1092: PetscObjectGetOptionsPrefix((PetscObject)mat,&mprefix);
1093: MatGetSubMatrices(mat,n_local,bjac->is,bjac->is,scall,&bjac->mat);
1094: }
1095: /* Return control to the user so that the submatrices can be modified (e.g., to apply
1096: different boundary conditions for the submatrices than for the global problem) */
1097: PCModifySubMatrices(pc,n_local,bjac->is,bjac->is,bjac->pmat,pc->modifysubmatricesP);
1099: PetscObjectGetOptionsPrefix((PetscObject)pmat,&pprefix);
1100: for (i=0; i<n_local; i++) {
1101: PetscLogObjectParent(pc,bjac->pmat[i]);
1102: PetscObjectSetOptionsPrefix((PetscObject)bjac->pmat[i],pprefix);
1103: if (pc->useAmat) {
1104: PetscLogObjectParent(pc,bjac->mat[i]);
1105: PetscObjectSetOptionsPrefix((PetscObject)bjac->mat[i],mprefix);
1106: KSPSetOperators(jac->ksp[i],bjac->mat[i],bjac->pmat[i],pc->flag);
1107: } else {
1108: KSPSetOperators(jac->ksp[i],bjac->pmat[i],bjac->pmat[i],pc->flag);
1109: }
1110: if (pc->setfromoptionscalled) {
1111: KSPSetFromOptions(jac->ksp[i]);
1112: }
1113: }
1114: return(0);
1115: }
1117: /* ---------------------------------------------------------------------------------------------*/
1118: /*
1119: These are for a single block with multiple processes;
1120: */
1123: static PetscErrorCode PCReset_BJacobi_Multiproc(PC pc)
1124: {
1125: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
1126: PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1127: PetscErrorCode ierr;
1130: VecDestroy(&mpjac->ysub);
1131: VecDestroy(&mpjac->xsub);
1132: MatDestroy(&mpjac->submats);
1133: if (jac->ksp) {KSPReset(jac->ksp[0]);}
1134: return(0);
1135: }
1139: static PetscErrorCode PCDestroy_BJacobi_Multiproc(PC pc)
1140: {
1141: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
1142: PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1143: PetscErrorCode ierr;
1146: PCReset_BJacobi_Multiproc(pc);
1147: KSPDestroy(&jac->ksp[0]);
1148: PetscFree(jac->ksp);
1149: PetscSubcommDestroy(&mpjac->psubcomm);
1151: PetscFree(mpjac);
1152: PetscFree(pc->data);
1153: return(0);
1154: }
1158: static PetscErrorCode PCApply_BJacobi_Multiproc(PC pc,Vec x,Vec y)
1159: {
1160: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
1161: PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1162: PetscErrorCode ierr;
1163: PetscScalar *xarray,*yarray;
1166: /* place x's and y's local arrays into xsub and ysub */
1167: VecGetArray(x,&xarray);
1168: VecGetArray(y,&yarray);
1169: VecPlaceArray(mpjac->xsub,xarray);
1170: VecPlaceArray(mpjac->ysub,yarray);
1172: /* apply preconditioner on each matrix block */
1173: PetscLogEventBegin(PC_ApplyOnMproc,jac->ksp[0],mpjac->xsub,mpjac->ysub,0);
1174: KSPSolve(jac->ksp[0],mpjac->xsub,mpjac->ysub);
1175: PetscLogEventEnd(PC_ApplyOnMproc,jac->ksp[0],mpjac->xsub,mpjac->ysub,0);
1177: VecResetArray(mpjac->xsub);
1178: VecResetArray(mpjac->ysub);
1179: VecRestoreArray(x,&xarray);
1180: VecRestoreArray(y,&yarray);
1181: return(0);
1182: }
1184: #include <petsc-private/matimpl.h>
1187: static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC pc)
1188: {
1189: PC_BJacobi *jac = (PC_BJacobi*)pc->data;
1190: PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data;
1191: PetscErrorCode ierr;
1192: PetscInt m,n;
1193: MPI_Comm comm,subcomm=0;
1194: const char *prefix;
1195: PetscBool wasSetup = PETSC_TRUE;
1198: PetscObjectGetComm((PetscObject)pc,&comm);
1199: if (jac->n_local > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only a single block in a subcommunicator is supported");
1200: jac->n_local = 1; /* currently only a single block is supported for a subcommunicator */
1201: if (!pc->setupcalled) {
1202: wasSetup = PETSC_FALSE;
1203: PetscNewLog(pc,PC_BJacobi_Multiproc,&mpjac);
1204: jac->data = (void*)mpjac;
1206: /* initialize datastructure mpjac */
1207: if (!jac->psubcomm) {
1208: /* Create default contiguous subcommunicatiors if user does not provide them */
1209: PetscSubcommCreate(comm,&jac->psubcomm);
1210: PetscSubcommSetNumber(jac->psubcomm,jac->n);
1211: PetscSubcommSetType(jac->psubcomm,PETSC_SUBCOMM_CONTIGUOUS);
1212: PetscLogObjectMemory(pc,sizeof(PetscSubcomm));
1213: }
1214: mpjac->psubcomm = jac->psubcomm;
1215: subcomm = mpjac->psubcomm->comm;
1217: /* Get matrix blocks of pmat */
1218: if (!pc->pmat->ops->getmultiprocblock) SETERRQ(PetscObjectComm((PetscObject)pc->pmat),PETSC_ERR_SUP,"No support for the requested operation");
1219: (*pc->pmat->ops->getmultiprocblock)(pc->pmat,subcomm,MAT_INITIAL_MATRIX,&mpjac->submats);
1221: /* create a new PC that processors in each subcomm have copy of */
1222: PetscMalloc(sizeof(KSP),&jac->ksp);
1223: KSPCreate(subcomm,&jac->ksp[0]);
1224: PetscObjectIncrementTabLevel((PetscObject)jac->ksp[0],(PetscObject)pc,1);
1225: PetscLogObjectParent(pc,jac->ksp[0]);
1226: KSPSetOperators(jac->ksp[0],mpjac->submats,mpjac->submats,pc->flag);
1227: KSPGetPC(jac->ksp[0],&mpjac->pc);
1229: PCGetOptionsPrefix(pc,&prefix);
1230: KSPSetOptionsPrefix(jac->ksp[0],prefix);
1231: KSPAppendOptionsPrefix(jac->ksp[0],"sub_");
1232: /*
1233: PetscMPIInt rank,subsize,subrank;
1234: MPI_Comm_rank(comm,&rank);
1235: MPI_Comm_size(subcomm,&subsize);
1236: MPI_Comm_rank(subcomm,&subrank);
1238: MatGetLocalSize(mpjac->submats,&m,NULL);
1239: MatGetSize(mpjac->submats,&n,NULL);
1240: PetscSynchronizedPrintf(comm,"[%d], sub-size %d,sub-rank %d\n",rank,subsize,subrank);
1241: PetscSynchronizedFlush(comm);
1242: */
1244: /* create dummy vectors xsub and ysub */
1245: MatGetLocalSize(mpjac->submats,&m,&n);
1246: VecCreateMPIWithArray(subcomm,1,n,PETSC_DECIDE,NULL,&mpjac->xsub);
1247: VecCreateMPIWithArray(subcomm,1,m,PETSC_DECIDE,NULL,&mpjac->ysub);
1248: PetscLogObjectParent(pc,mpjac->xsub);
1249: PetscLogObjectParent(pc,mpjac->ysub);
1251: pc->ops->reset = PCReset_BJacobi_Multiproc;
1252: pc->ops->destroy = PCDestroy_BJacobi_Multiproc;
1253: pc->ops->apply = PCApply_BJacobi_Multiproc;
1254: } else { /* pc->setupcalled */
1255: subcomm = mpjac->psubcomm->comm;
1256: if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
1257: /* destroy old matrix blocks, then get new matrix blocks */
1258: if (mpjac->submats) {MatDestroy(&mpjac->submats);}
1259: (*pc->pmat->ops->getmultiprocblock)(pc->pmat,subcomm,MAT_INITIAL_MATRIX,&mpjac->submats);
1260: } else {
1261: (*pc->pmat->ops->getmultiprocblock)(pc->pmat,subcomm,MAT_REUSE_MATRIX,&mpjac->submats);
1262: }
1263: KSPSetOperators(jac->ksp[0],mpjac->submats,mpjac->submats,pc->flag);
1264: }
1266: if (!wasSetup && pc->setfromoptionscalled) {
1267: KSPSetFromOptions(jac->ksp[0]);
1268: }
1269: return(0);
1270: }