Actual source code: redundant.c
1: #define PETSCKSP_DLL
3: /*
4: This file defines a "solve the problem redundantly on each subgroup of processor" preconditioner.
5: */
6: #include private/pcimpl.h
7: #include petscksp.h
9: typedef struct {
10: KSP ksp;
11: PC pc; /* actual preconditioner used on each processor */
12: Vec xsub,ysub; /* vectors of a subcommunicator to hold parallel vectors of ((PetscObject)pc)->comm */
13: Vec xdup,ydup; /* parallel vector that congregates xsub or ysub facilitating vector scattering */
14: Mat pmats; /* matrix and optional preconditioner matrix belong to a subcommunicator */
15: VecScatter scatterin,scatterout; /* scatter used to move all values to each processor group (subcommunicator) */
16: PetscTruth useparallelmat;
17: PetscSubcomm psubcomm;
18: PetscInt nsubcomm; /* num of data structure PetscSubcomm */
19: } PC_Redundant;
23: static PetscErrorCode PCView_Redundant(PC pc,PetscViewer viewer)
24: {
25: PC_Redundant *red = (PC_Redundant*)pc->data;
27: PetscMPIInt rank;
28: PetscTruth iascii,isstring;
29: PetscViewer subviewer;
32: MPI_Comm_rank(((PetscObject)pc)->comm,&rank);
33: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
34: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
35: if (iascii) {
36: if (!red->psubcomm) {
37: PetscViewerASCIIPrintf(viewer," Redundant preconditioner: Not yet setup\n");
38: } else {
39: PetscViewerASCIIPrintf(viewer," Redundant preconditioner: First (color=0) of %D PCs follows\n",red->nsubcomm);
40: PetscViewerGetSubcomm(viewer,((PetscObject)red->pc)->comm,&subviewer);
41: if (red->psubcomm->color) { /* only view first redundant pc */
42: PetscViewerASCIIPushTab(viewer);
43: KSPView(red->ksp,subviewer);
44: PetscViewerASCIIPopTab(viewer);
45: }
46: PetscViewerRestoreSubcomm(viewer,((PetscObject)red->pc)->comm,&subviewer);
47: }
48: } else if (isstring) {
49: PetscViewerStringSPrintf(viewer," Redundant solver preconditioner");
50: } else {
51: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PC redundant",((PetscObject)viewer)->type_name);
52: }
53: return(0);
54: }
58: static PetscErrorCode PCSetUp_Redundant(PC pc)
59: {
60: PC_Redundant *red = (PC_Redundant*)pc->data;
62: PetscInt mstart,mend,mlocal,m,mlocal_sub,rstart_sub,rend_sub,mloc_sub;
63: PetscMPIInt size;
64: MatReuse reuse = MAT_INITIAL_MATRIX;
65: MatStructure str = DIFFERENT_NONZERO_PATTERN;
66: MPI_Comm comm = ((PetscObject)pc)->comm,subcomm;
67: Vec vec;
68: PetscMPIInt subsize,subrank;
69: const char *prefix;
70: const PetscInt *range;
73: MatGetVecs(pc->pmat,&vec,0);
74: VecGetSize(vec,&m);
76: if (!pc->setupcalled) {
77: if (!red->psubcomm) {
78: PetscSubcommCreate(comm,red->nsubcomm,&red->psubcomm);
79: PetscLogObjectMemory(pc,sizeof(PetscSubcomm));
81: /* create a new PC that processors in each subcomm have copy of */
82: subcomm = red->psubcomm->comm;
83: KSPCreate(subcomm,&red->ksp);
84: PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);
85: PetscLogObjectParent(pc,red->ksp);
86: KSPSetType(red->ksp,KSPPREONLY);
87: KSPGetPC(red->ksp,&red->pc);
88: PCSetType(red->pc,PCLU);
90: PCGetOptionsPrefix(pc,&prefix);
91: KSPSetOptionsPrefix(red->ksp,prefix);
92: KSPAppendOptionsPrefix(red->ksp,"redundant_");
93: } else {
94: subcomm = red->psubcomm->comm;
95: }
97: /* create working vectors xsub/ysub and xdup/ydup */
98: VecGetLocalSize(vec,&mlocal);
99: VecGetOwnershipRange(vec,&mstart,&mend);
101: /* get local size of xsub/ysub */
102: MPI_Comm_size(subcomm,&subsize);
103: MPI_Comm_rank(subcomm,&subrank);
104: MatGetOwnershipRanges(pc->pmat,&range);
105: rstart_sub = range[red->psubcomm->n*subrank]; /* rstart in xsub/ysub */
106: if (subrank+1 < subsize){
107: rend_sub = range[red->psubcomm->n*(subrank+1)];
108: } else {
109: rend_sub = m;
110: }
111: mloc_sub = rend_sub - rstart_sub;
112: VecCreateMPI(subcomm,mloc_sub,PETSC_DECIDE,&red->ysub);
113: /* create xsub with empty local arrays, because xdup's arrays will be placed into it */
114: VecCreateMPIWithArray(subcomm,mloc_sub,PETSC_DECIDE,PETSC_NULL,&red->xsub);
116: /* create xdup and ydup. ydup has empty local arrays because ysub's arrays will be place into it.
117: Note: we use communicator dupcomm, not ((PetscObject)pc)->comm! */
118: VecCreateMPI(red->psubcomm->dupparent,mloc_sub,PETSC_DECIDE,&red->xdup);
119: VecCreateMPIWithArray(red->psubcomm->dupparent,mloc_sub,PETSC_DECIDE,PETSC_NULL,&red->ydup);
120:
121: /* create vec scatters */
122: if (!red->scatterin){
123: IS is1,is2;
124: PetscInt *idx1,*idx2,i,j,k;
126: PetscMalloc2(red->psubcomm->n*mlocal,PetscInt,&idx1,red->psubcomm->n*mlocal,PetscInt,&idx2);
127: j = 0;
128: for (k=0; k<red->psubcomm->n; k++){
129: for (i=mstart; i<mend; i++){
130: idx1[j] = i;
131: idx2[j++] = i + m*k;
132: }
133: }
134: ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx1,&is1);
135: ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx2,&is2);
136: VecScatterCreate(vec,is1,red->xdup,is2,&red->scatterin);
137: ISDestroy(is1);
138: ISDestroy(is2);
140: ISCreateStride(comm,mlocal,mstart+ red->psubcomm->color*m,1,&is1);
141: ISCreateStride(comm,mlocal,mstart,1,&is2);
142: VecScatterCreate(red->xdup,is1,vec,is2,&red->scatterout);
143: ISDestroy(is1);
144: ISDestroy(is2);
145: PetscFree2(idx1,idx2);
146: }
147: }
148: VecDestroy(vec);
150: /* if pmatrix set by user is sequential then we do not need to gather the parallel matrix */
151: MPI_Comm_size(comm,&size);
152: if (size == 1) {
153: red->useparallelmat = PETSC_FALSE;
154: }
156: if (red->useparallelmat) {
157: if (pc->setupcalled == 1 && pc->flag == DIFFERENT_NONZERO_PATTERN) {
158: /* destroy old matrices */
159: if (red->pmats) {
160: MatDestroy(red->pmats);
161: }
162: } else if (pc->setupcalled == 1) {
163: reuse = MAT_REUSE_MATRIX;
164: str = SAME_NONZERO_PATTERN;
165: }
166:
167: /* grab the parallel matrix and put it into processors of a subcomminicator */
168: /*--------------------------------------------------------------------------*/
169: VecGetLocalSize(red->ysub,&mlocal_sub);
170: MatGetRedundantMatrix(pc->pmat,red->psubcomm->n,red->psubcomm->comm,mlocal_sub,reuse,&red->pmats);
171: /* tell PC of the subcommunicator its operators */
172: KSPSetOperators(red->ksp,red->pmats,red->pmats,str);
173: } else {
174: KSPSetOperators(red->ksp,pc->mat,pc->pmat,pc->flag);
175: }
176: if (pc->setfromoptionscalled){
177: KSPSetFromOptions(red->ksp);
178: }
179: KSPSetUp(red->ksp);
180: return(0);
181: }
185: static PetscErrorCode PCApply_Redundant(PC pc,Vec x,Vec y)
186: {
187: PC_Redundant *red = (PC_Redundant*)pc->data;
189: PetscScalar *array;
192: /* scatter x to xdup */
193: VecScatterBegin(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);
194: VecScatterEnd(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);
195:
196: /* place xdup's local array into xsub */
197: VecGetArray(red->xdup,&array);
198: VecPlaceArray(red->xsub,(const PetscScalar*)array);
200: /* apply preconditioner on each processor */
201: PCApply(red->pc,red->xsub,red->ysub);
202: VecResetArray(red->xsub);
203: VecRestoreArray(red->xdup,&array);
204:
205: /* place ysub's local array into ydup */
206: VecGetArray(red->ysub,&array);
207: VecPlaceArray(red->ydup,(const PetscScalar*)array);
209: /* scatter ydup to y */
210: VecScatterBegin(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);
211: VecScatterEnd(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);
212: VecResetArray(red->ydup);
213: VecRestoreArray(red->ysub,&array);
214: return(0);
215: }
219: static PetscErrorCode PCDestroy_Redundant(PC pc)
220: {
221: PC_Redundant *red = (PC_Redundant*)pc->data;
225: if (red->scatterin) {VecScatterDestroy(red->scatterin);}
226: if (red->scatterout) {VecScatterDestroy(red->scatterout);}
227: if (red->ysub) {VecDestroy(red->ysub);}
228: if (red->xsub) {VecDestroy(red->xsub);}
229: if (red->xdup) {VecDestroy(red->xdup);}
230: if (red->ydup) {VecDestroy(red->ydup);}
231: if (red->pmats) {
232: MatDestroy(red->pmats);
233: }
234: if (red->psubcomm) {PetscSubcommDestroy(red->psubcomm);}
235: if (red->ksp) {KSPDestroy(red->ksp);}
236: PetscFree(red);
237: return(0);
238: }
242: static PetscErrorCode PCSetFromOptions_Redundant(PC pc)
243: {
245: PC_Redundant *red = (PC_Redundant*)pc->data;
248: PetscOptionsHead("Redundant options");
249: PetscOptionsInt("-pc_redundant_number","Number of redundant pc","PCRedundantSetNumber",red->nsubcomm,&red->nsubcomm,0);
250: PetscOptionsTail();
251: return(0);
252: }
257: PetscErrorCode PCRedundantSetNumber_Redundant(PC pc,PetscInt nreds)
258: {
259: PC_Redundant *red = (PC_Redundant*)pc->data;
262: red->nsubcomm = nreds;
263: return(0);
264: }
269: /*@
270: PCRedundantSetNumber - Sets the number of redundant preconditioner contexts.
272: Collective on PC
274: Input Parameters:
275: + pc - the preconditioner context
276: - nredundant - number of redundant preconditioner contexts; for example if you are using 64 MPI processes and
277: use an nredundant of 4 there will be 4 parallel solves each on 16 = 64/4 processes.
279: Level: advanced
281: .keywords: PC, redundant solve
282: @*/
283: PetscErrorCode PCRedundantSetNumber(PC pc,PetscInt nredundant)
284: {
285: PetscErrorCode ierr,(*f)(PC,PetscInt);
289: if (nredundant <= 0) SETERRQ1(PETSC_ERR_ARG_WRONG, "num of redundant pc %D must be positive",nredundant);
290: PetscObjectQueryFunction((PetscObject)pc,"PCRedundantSetNumber_C",(void (**)(void))&f);
291: if (f) {
292: (*f)(pc,nredundant);
293: }
294: return(0);
295: }
300: PetscErrorCode PCRedundantSetScatter_Redundant(PC pc,VecScatter in,VecScatter out)
301: {
302: PC_Redundant *red = (PC_Redundant*)pc->data;
306: PetscObjectReference((PetscObject)in);
307: if (red->scatterin) { VecScatterDestroy(red->scatterin); }
308: red->scatterin = in;
309: PetscObjectReference((PetscObject)out);
310: if (red->scatterout) { VecScatterDestroy(red->scatterout); }
311: red->scatterout = out;
312: return(0);
313: }
318: /*@
319: PCRedundantSetScatter - Sets the scatter used to copy values into the
320: redundant local solve and the scatter to move them back into the global
321: vector.
323: Collective on PC
325: Input Parameters:
326: + pc - the preconditioner context
327: . in - the scatter to move the values in
328: - out - the scatter to move them out
330: Level: advanced
332: .keywords: PC, redundant solve
333: @*/
334: PetscErrorCode PCRedundantSetScatter(PC pc,VecScatter in,VecScatter out)
335: {
336: PetscErrorCode ierr,(*f)(PC,VecScatter,VecScatter);
342: PetscObjectQueryFunction((PetscObject)pc,"PCRedundantSetScatter_C",(void (**)(void))&f);
343: if (f) {
344: (*f)(pc,in,out);
345: }
346: return(0);
347: }
352: PetscErrorCode PCRedundantGetPC_Redundant(PC pc,PC *innerpc)
353: {
355: PC_Redundant *red = (PC_Redundant*)pc->data;
356: MPI_Comm comm,subcomm;
357: const char *prefix;
360: if (!red->psubcomm) {
361: PetscObjectGetComm((PetscObject)pc,&comm);
362: PetscSubcommCreate(comm,red->nsubcomm,&red->psubcomm);
363: PetscLogObjectMemory(pc,sizeof(PetscSubcomm));
365: /* create a new PC that processors in each subcomm have copy of */
366: subcomm = red->psubcomm->comm;
367: KSPCreate(subcomm,&red->ksp);
368: PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);
369: PetscLogObjectParent(pc,red->ksp);
370: KSPSetType(red->ksp,KSPPREONLY);
371: KSPGetPC(red->ksp,&red->pc);
372: PCSetType(red->pc,PCLU);
374: PCGetOptionsPrefix(pc,&prefix);
375: KSPSetOptionsPrefix(red->ksp,prefix);
376: KSPAppendOptionsPrefix(red->ksp,"redundant_");
377: }
379: KSPGetPC(red->ksp,innerpc);
380: return(0);
381: }
386: /*@
387: PCRedundantGetPC - Gets the sequential PC created by the redundant PC.
389: Not Collective
391: Input Parameter:
392: . pc - the preconditioner context
394: Output Parameter:
395: . innerpc - the sequential PC
397: Level: advanced
399: .keywords: PC, redundant solve
400: @*/
401: PetscErrorCode PCRedundantGetPC(PC pc,PC *innerpc)
402: {
403: PetscErrorCode ierr,(*f)(PC,PC*);
408: PetscObjectQueryFunction((PetscObject)pc,"PCRedundantGetPC_C",(void (**)(void))&f);
409: if (f) {
410: (*f)(pc,innerpc);
411: }
412: return(0);
413: }
418: PetscErrorCode PCRedundantGetOperators_Redundant(PC pc,Mat *mat,Mat *pmat)
419: {
420: PC_Redundant *red = (PC_Redundant*)pc->data;
423: if (mat) *mat = red->pmats;
424: if (pmat) *pmat = red->pmats;
425: return(0);
426: }
431: /*@
432: PCRedundantGetOperators - gets the sequential matrix and preconditioner matrix
434: Not Collective
436: Input Parameter:
437: . pc - the preconditioner context
439: Output Parameters:
440: + mat - the matrix
441: - pmat - the (possibly different) preconditioner matrix
443: Level: advanced
445: .keywords: PC, redundant solve
446: @*/
447: PetscErrorCode PCRedundantGetOperators(PC pc,Mat *mat,Mat *pmat)
448: {
449: PetscErrorCode ierr,(*f)(PC,Mat*,Mat*);
455: PetscObjectQueryFunction((PetscObject)pc,"PCRedundantGetOperators_C",(void (**)(void))&f);
456: if (f) {
457: (*f)(pc,mat,pmat);
458: }
459: return(0);
460: }
462: /* -------------------------------------------------------------------------------------*/
463: /*MC
464: PCREDUNDANT - Runs a preconditioner for the entire problem on subgroups of processors
466: Options for the redundant preconditioners can be set with -redundant_pc_xxx
468: Options Database:
469: . -pc_redundant_number <n> - number of redundant solves, for example if you are using 64 MPI processes and
470: use an n of 4 there will be 4 parallel solves each on 16 = 64/4 processes.
472: Level: intermediate
474: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCRedundantSetScatter(),
475: PCRedundantGetPC(), PCRedundantGetOperators(), PCRedundantSetNumber()
476: M*/
481: PetscErrorCode PCCreate_Redundant(PC pc)
482: {
484: PC_Redundant *red;
485: PetscMPIInt size;
486:
488: PetscNewLog(pc,PC_Redundant,&red);
489: MPI_Comm_size(((PetscObject)pc)->comm,&size);
490: red->nsubcomm = size;
491: red->useparallelmat = PETSC_TRUE;
492: pc->data = (void*)red;
494: pc->ops->apply = PCApply_Redundant;
495: pc->ops->applytranspose = 0;
496: pc->ops->setup = PCSetUp_Redundant;
497: pc->ops->destroy = PCDestroy_Redundant;
498: pc->ops->setfromoptions = PCSetFromOptions_Redundant;
499: pc->ops->view = PCView_Redundant;
500: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantSetScatter_C","PCRedundantSetScatter_Redundant",
501: PCRedundantSetScatter_Redundant);
502: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantSetNumber_C","PCRedundantSetNumber_Redundant",
503: PCRedundantSetNumber_Redundant);
504: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantGetPC_C","PCRedundantGetPC_Redundant",
505: PCRedundantGetPC_Redundant);
506: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantGetOperators_C","PCRedundantGetOperators_Redundant",
507: PCRedundantGetOperators_Redundant);
508: return(0);
509: }