Actual source code: redundant.c
petsc-dev 2014-02-02
2: /*
3: This file defines a "solve the problem redundantly on each subgroup of processor" preconditioner.
4: */
5: #include <petsc-private/pcimpl.h>
6: #include <petscksp.h> /*I "petscksp.h" I*/
8: typedef struct {
9: KSP ksp;
10: PC pc; /* actual preconditioner used on each processor */
11: Vec xsub,ysub; /* vectors of a subcommunicator to hold parallel vectors of PetscObjectComm((PetscObject)pc) */
12: Vec xdup,ydup; /* parallel vector that congregates xsub or ysub facilitating vector scattering */
13: Mat pmats; /* matrix and optional preconditioner matrix belong to a subcommunicator */
14: VecScatter scatterin,scatterout; /* scatter used to move all values to each processor group (subcommunicator) */
15: PetscBool useparallelmat;
16: PetscSubcomm psubcomm;
17: PetscInt nsubcomm; /* num of data structure PetscSubcomm */
18: } PC_Redundant;
22: static PetscErrorCode PCView_Redundant(PC pc,PetscViewer viewer)
23: {
24: PC_Redundant *red = (PC_Redundant*)pc->data;
26: PetscBool iascii,isstring;
27: PetscViewer subviewer;
30: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
31: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
32: if (iascii) {
33: if (!red->psubcomm) {
34: PetscViewerASCIIPrintf(viewer," Redundant preconditioner: Not yet setup\n");
35: } else {
36: PetscViewerASCIIPrintf(viewer," Redundant preconditioner: First (color=0) of %D PCs follows\n",red->nsubcomm);
37: PetscViewerGetSubcomm(viewer,((PetscObject)red->pc)->comm,&subviewer);
38: if (!red->psubcomm->color) { /* only view first redundant pc */
39: PetscViewerASCIIPushTab(viewer);
40: KSPView(red->ksp,subviewer);
41: PetscViewerASCIIPopTab(viewer);
42: }
43: PetscViewerRestoreSubcomm(viewer,((PetscObject)red->pc)->comm,&subviewer);
44: }
45: } else if (isstring) {
46: PetscViewerStringSPrintf(viewer," Redundant solver preconditioner");
47: }
48: return(0);
49: }
53: static PetscErrorCode PCSetUp_Redundant(PC pc)
54: {
55: PC_Redundant *red = (PC_Redundant*)pc->data;
57: PetscInt mstart,mend,mlocal,M;
58: PetscMPIInt size;
59: MPI_Comm comm,subcomm;
60: Vec x;
61: const char *prefix;
64: PetscObjectGetComm((PetscObject)pc,&comm);
65:
66: /* if pmatrix set by user is sequential then we do not need to gather the parallel matrix */
67: MPI_Comm_size(comm,&size);
68: if (size == 1) red->useparallelmat = PETSC_FALSE;
70: if (!pc->setupcalled) {
71: PetscInt mloc_sub;
72: if (!red->psubcomm) {
73: PetscSubcommCreate(comm,&red->psubcomm);
74: PetscSubcommSetNumber(red->psubcomm,red->nsubcomm);
75: PetscSubcommSetType(red->psubcomm,PETSC_SUBCOMM_CONTIGUOUS);
76: /* enable runtime switch of psubcomm type, e.g., '-psubcomm_type interlaced */
77: PetscSubcommSetFromOptions(red->psubcomm);
78: PetscLogObjectMemory((PetscObject)pc,sizeof(PetscSubcomm));
80: /* create a new PC that processors in each subcomm have copy of */
81: subcomm = red->psubcomm->comm;
83: KSPCreate(subcomm,&red->ksp);
84: PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);
85: PetscLogObjectParent((PetscObject)pc,(PetscObject)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: if (red->useparallelmat) {
98: /* grab the parallel matrix and put it into processors of a subcomminicator */
99: MatGetRedundantMatrix(pc->pmat,red->psubcomm->n,subcomm,MAT_INITIAL_MATRIX,&red->pmats);
100: KSPSetOperators(red->ksp,red->pmats,red->pmats,DIFFERENT_NONZERO_PATTERN);
101:
102: /* get working vectors xsub and ysub */
103: MatGetVecs(red->pmats,&red->xsub,&red->ysub);
105: /* create working vectors xdup and ydup.
106: xdup concatenates all xsub's contigously to form a mpi vector over dupcomm (see PetscSubcommCreate_interlaced())
107: ydup concatenates all ysub and has empty local arrays because ysub's arrays will be place into it.
108: Note: we use communicator dupcomm, not PetscObjectComm((PetscObject)pc)! */
109: MatGetLocalSize(red->pmats,&mloc_sub,NULL);
110: VecCreateMPI(red->psubcomm->dupparent,mloc_sub,PETSC_DECIDE,&red->xdup);
111: VecCreateMPIWithArray(red->psubcomm->dupparent,1,mloc_sub,PETSC_DECIDE,NULL,&red->ydup);
113: /* create vecscatters */
114: if (!red->scatterin) { /* efficiency of scatterin is independent from psubcomm_type! */
115: IS is1,is2;
116: PetscInt *idx1,*idx2,i,j,k;
118: MatGetVecs(pc->pmat,&x,0);
119: VecGetSize(x,&M);
120: VecGetOwnershipRange(x,&mstart,&mend);
121: mlocal = mend - mstart;
122: PetscMalloc2(red->psubcomm->n*mlocal,&idx1,red->psubcomm->n*mlocal,&idx2);
123: j = 0;
124: for (k=0; k<red->psubcomm->n; k++) {
125: for (i=mstart; i<mend; i++) {
126: idx1[j] = i;
127: idx2[j++] = i + M*k;
128: }
129: }
130: ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx1,PETSC_COPY_VALUES,&is1);
131: ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx2,PETSC_COPY_VALUES,&is2);
132: VecScatterCreate(x,is1,red->xdup,is2,&red->scatterin);
133: ISDestroy(&is1);
134: ISDestroy(&is2);
136: /* Impl below is good for PETSC_SUBCOMM_INTERLACED (no inter-process communication) and PETSC_SUBCOMM_CONTIGUOUS (communication within subcomm) */
137: ISCreateStride(comm,mlocal,mstart+ red->psubcomm->color*M,1,&is1);
138: ISCreateStride(comm,mlocal,mstart,1,&is2);
139: VecScatterCreate(red->xdup,is1,x,is2,&red->scatterout);
140: ISDestroy(&is1);
141: ISDestroy(&is2);
142: PetscFree2(idx1,idx2);
143: VecDestroy(&x);
144: }
145: } else { /* !red->useparallelmat */
146: KSPSetOperators(red->ksp,pc->mat,pc->pmat,pc->flag);
147: }
148: } else { /* pc->setupcalled */
149: if (red->useparallelmat) {
150: MatReuse reuse;
151: /* grab the parallel matrix and put it into processors of a subcomminicator */
152: /*--------------------------------------------------------------------------*/
153: if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
154: /* destroy old matrices */
155: MatDestroy(&red->pmats);
156: reuse = MAT_INITIAL_MATRIX;
157: } else {
158: reuse = MAT_REUSE_MATRIX;
159: }
160: MatGetRedundantMatrix(pc->pmat,red->psubcomm->n,red->psubcomm->comm,reuse,&red->pmats);
161: KSPSetOperators(red->ksp,red->pmats,red->pmats,pc->flag);
162: } else { /* !red->useparallelmat */
163: KSPSetOperators(red->ksp,pc->mat,pc->pmat,pc->flag);
164: }
165: }
167: if (pc->setfromoptionscalled) {
168: KSPSetFromOptions(red->ksp);
169: }
170: KSPSetUp(red->ksp);
171: return(0);
172: }
176: static PetscErrorCode PCApply_Redundant(PC pc,Vec x,Vec y)
177: {
178: PC_Redundant *red = (PC_Redundant*)pc->data;
180: PetscScalar *array;
183: if (!red->useparallelmat) {
184: KSPSolve(red->ksp,x,y);
185: return(0);
186: }
188: /* scatter x to xdup */
189: VecScatterBegin(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);
190: VecScatterEnd(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);
192: /* place xdup's local array into xsub */
193: VecGetArray(red->xdup,&array);
194: VecPlaceArray(red->xsub,(const PetscScalar*)array);
196: /* apply preconditioner on each processor */
197: KSPSolve(red->ksp,red->xsub,red->ysub);
198: VecResetArray(red->xsub);
199: VecRestoreArray(red->xdup,&array);
201: /* place ysub's local array into ydup */
202: VecGetArray(red->ysub,&array);
203: VecPlaceArray(red->ydup,(const PetscScalar*)array);
205: /* scatter ydup to y */
206: VecScatterBegin(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);
207: VecScatterEnd(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);
208: VecResetArray(red->ydup);
209: VecRestoreArray(red->ysub,&array);
210: return(0);
211: }
215: static PetscErrorCode PCReset_Redundant(PC pc)
216: {
217: PC_Redundant *red = (PC_Redundant*)pc->data;
221: if (red->useparallelmat) {
222: VecScatterDestroy(&red->scatterin);
223: VecScatterDestroy(&red->scatterout);
224: VecDestroy(&red->ysub);
225: VecDestroy(&red->xsub);
226: VecDestroy(&red->xdup);
227: VecDestroy(&red->ydup);
228: }
229: MatDestroy(&red->pmats);
230: KSPReset(red->ksp);
231: return(0);
232: }
236: static PetscErrorCode PCDestroy_Redundant(PC pc)
237: {
238: PC_Redundant *red = (PC_Redundant*)pc->data;
242: PCReset_Redundant(pc);
243: KSPDestroy(&red->ksp);
244: PetscSubcommDestroy(&red->psubcomm);
245: PetscFree(pc->data);
246: return(0);
247: }
251: static PetscErrorCode PCSetFromOptions_Redundant(PC pc)
252: {
254: PC_Redundant *red = (PC_Redundant*)pc->data;
257: PetscOptionsHead("Redundant options");
258: PetscOptionsInt("-pc_redundant_number","Number of redundant pc","PCRedundantSetNumber",red->nsubcomm,&red->nsubcomm,0);
259: PetscOptionsTail();
260: return(0);
261: }
265: static PetscErrorCode PCRedundantSetNumber_Redundant(PC pc,PetscInt nreds)
266: {
267: PC_Redundant *red = (PC_Redundant*)pc->data;
270: red->nsubcomm = nreds;
271: return(0);
272: }
276: /*@
277: PCRedundantSetNumber - Sets the number of redundant preconditioner contexts.
279: Logically Collective on PC
281: Input Parameters:
282: + pc - the preconditioner context
283: - nredundant - number of redundant preconditioner contexts; for example if you are using 64 MPI processes and
284: use an nredundant of 4 there will be 4 parallel solves each on 16 = 64/4 processes.
286: Level: advanced
288: .keywords: PC, redundant solve
289: @*/
290: PetscErrorCode PCRedundantSetNumber(PC pc,PetscInt nredundant)
291: {
296: if (nredundant <= 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "num of redundant pc %D must be positive",nredundant);
297: PetscTryMethod(pc,"PCRedundantSetNumber_C",(PC,PetscInt),(pc,nredundant));
298: return(0);
299: }
303: static PetscErrorCode PCRedundantSetScatter_Redundant(PC pc,VecScatter in,VecScatter out)
304: {
305: PC_Redundant *red = (PC_Redundant*)pc->data;
309: PetscObjectReference((PetscObject)in);
310: VecScatterDestroy(&red->scatterin);
312: red->scatterin = in;
314: PetscObjectReference((PetscObject)out);
315: VecScatterDestroy(&red->scatterout);
316: red->scatterout = out;
317: return(0);
318: }
322: /*@
323: PCRedundantSetScatter - Sets the scatter used to copy values into the
324: redundant local solve and the scatter to move them back into the global
325: vector.
327: Logically Collective on PC
329: Input Parameters:
330: + pc - the preconditioner context
331: . in - the scatter to move the values in
332: - out - the scatter to move them out
334: Level: advanced
336: .keywords: PC, redundant solve
337: @*/
338: PetscErrorCode PCRedundantSetScatter(PC pc,VecScatter in,VecScatter out)
339: {
346: PetscTryMethod(pc,"PCRedundantSetScatter_C",(PC,VecScatter,VecScatter),(pc,in,out));
347: return(0);
348: }
352: static PetscErrorCode PCRedundantGetKSP_Redundant(PC pc,KSP *innerksp)
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->psubcomm);
363: PetscSubcommSetNumber(red->psubcomm,red->nsubcomm);
364: PetscSubcommSetType(red->psubcomm,PETSC_SUBCOMM_INTERLACED);
365: PetscLogObjectMemory((PetscObject)pc,sizeof(PetscSubcomm));
367: /* create a new PC that processors in each subcomm have copy of */
368: subcomm = red->psubcomm->comm;
370: KSPCreate(subcomm,&red->ksp);
371: PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);
372: PetscLogObjectParent((PetscObject)pc,(PetscObject)red->ksp);
373: KSPSetType(red->ksp,KSPPREONLY);
374: KSPGetPC(red->ksp,&red->pc);
375: PCSetType(red->pc,PCLU);
377: PCGetOptionsPrefix(pc,&prefix);
378: KSPSetOptionsPrefix(red->ksp,prefix);
379: KSPAppendOptionsPrefix(red->ksp,"redundant_");
380: }
381: *innerksp = red->ksp;
382: return(0);
383: }
387: /*@
388: PCRedundantGetKSP - Gets the less parallel KSP created by the redundant PC.
390: Not Collective
392: Input Parameter:
393: . pc - the preconditioner context
395: Output Parameter:
396: . innerksp - the KSP on the smaller set of processes
398: Level: advanced
400: .keywords: PC, redundant solve
401: @*/
402: PetscErrorCode PCRedundantGetKSP(PC pc,KSP *innerksp)
403: {
409: PetscTryMethod(pc,"PCRedundantGetKSP_C",(PC,KSP*),(pc,innerksp));
410: return(0);
411: }
415: static PetscErrorCode PCRedundantGetOperators_Redundant(PC pc,Mat *mat,Mat *pmat)
416: {
417: PC_Redundant *red = (PC_Redundant*)pc->data;
420: if (mat) *mat = red->pmats;
421: if (pmat) *pmat = red->pmats;
422: return(0);
423: }
427: /*@
428: PCRedundantGetOperators - gets the sequential matrix and preconditioner matrix
430: Not Collective
432: Input Parameter:
433: . pc - the preconditioner context
435: Output Parameters:
436: + mat - the matrix
437: - pmat - the (possibly different) preconditioner matrix
439: Level: advanced
441: .keywords: PC, redundant solve
442: @*/
443: PetscErrorCode PCRedundantGetOperators(PC pc,Mat *mat,Mat *pmat)
444: {
451: PetscTryMethod(pc,"PCRedundantGetOperators_C",(PC,Mat*,Mat*),(pc,mat,pmat));
452: return(0);
453: }
455: /* -------------------------------------------------------------------------------------*/
456: /*MC
457: PCREDUNDANT - Runs a KSP solver with preconditioner for the entire problem on subgroups of processors
459: Options for the redundant preconditioners can be set with -redundant_pc_xxx for the redundant KSP with -redundant_ksp_xxx
461: Options Database:
462: . -pc_redundant_number <n> - number of redundant solves, for example if you are using 64 MPI processes and
463: use an n of 4 there will be 4 parallel solves each on 16 = 64/4 processes.
465: Level: intermediate
467: Notes: The default KSP is preonly and the default PC is LU.
469: Developer Notes: Note that PCSetInitialGuessNonzero() is not used by this class but likely should be.
471: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCRedundantSetScatter(),
472: PCRedundantGetKSP(), PCRedundantGetOperators(), PCRedundantSetNumber()
473: M*/
477: PETSC_EXTERN PetscErrorCode PCCreate_Redundant(PC pc)
478: {
480: PC_Redundant *red;
481: PetscMPIInt size;
484: PetscNewLog(pc,&red);
485: MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
487: red->nsubcomm = size;
488: red->useparallelmat = PETSC_TRUE;
489: pc->data = (void*)red;
491: pc->ops->apply = PCApply_Redundant;
492: pc->ops->applytranspose = 0;
493: pc->ops->setup = PCSetUp_Redundant;
494: pc->ops->destroy = PCDestroy_Redundant;
495: pc->ops->reset = PCReset_Redundant;
496: pc->ops->setfromoptions = PCSetFromOptions_Redundant;
497: pc->ops->view = PCView_Redundant;
499: PetscObjectComposeFunction((PetscObject)pc,"PCRedundantSetScatter_C",PCRedundantSetScatter_Redundant);
500: PetscObjectComposeFunction((PetscObject)pc,"PCRedundantSetNumber_C",PCRedundantSetNumber_Redundant);
501: PetscObjectComposeFunction((PetscObject)pc,"PCRedundantGetKSP_C",PCRedundantGetKSP_Redundant);
502: PetscObjectComposeFunction((PetscObject)pc,"PCRedundantGetOperators_C",PCRedundantGetOperators_Redundant);
503: return(0);
504: }