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 sviewer,subviewer;
30: PetscInt color = red->psubcomm->color;
33: MPI_Comm_rank(((PetscObject)pc)->comm,&rank);
34: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
35: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
36: if (iascii) {
37: PetscViewerASCIIPrintf(viewer," Redundant preconditioner: First (color=0) of %D PCs follows\n",red->nsubcomm);
38: PetscViewerGetSubcomm(viewer,((PetscObject)red->pc)->comm,&subviewer);
39: if (!color) { /* only view first redundant pc */
40: PetscViewerASCIIPushTab(viewer);
41: KSPView(red->ksp,subviewer);
42: PetscViewerASCIIPopTab(viewer);
43: }
44: PetscViewerRestoreSubcomm(viewer,((PetscObject)red->pc)->comm,&subviewer);
45: } else if (isstring) { /* not test it yet! */
46: PetscViewerStringSPrintf(viewer," Redundant solver preconditioner");
47: PetscViewerGetSingleton(viewer,&sviewer);
48: if (!rank) {
49: KSPView(red->ksp,sviewer);
50: }
51: PetscViewerRestoreSingleton(viewer,&sviewer);
52: } else {
53: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PC redundant",((PetscObject)viewer)->type_name);
54: }
55: return(0);
56: }
58: #include private/matimpl.h
61: static PetscErrorCode PCSetUp_Redundant(PC pc)
62: {
63: PC_Redundant *red = (PC_Redundant*)pc->data;
65: PetscInt mstart,mend,mlocal,m,mlocal_sub,rstart_sub,rend_sub,mloc_sub;
66: PetscMPIInt size;
67: MatReuse reuse = MAT_INITIAL_MATRIX;
68: MatStructure str = DIFFERENT_NONZERO_PATTERN;
69: MPI_Comm comm = ((PetscObject)pc)->comm,subcomm;
70: Vec vec;
71: PetscMPIInt subsize,subrank;
72: const char *prefix;
75: MatGetVecs(pc->pmat,&vec,0);
76: VecGetSize(vec,&m);
78: if (!pc->setupcalled) {
79: if (!red->psubcomm) {
80: PetscSubcommCreate(comm,red->nsubcomm,&red->psubcomm);
81: PetscLogObjectMemory(pc,sizeof(PetscSubcomm));
83: /* create a new PC that processors in each subcomm have copy of */
84: subcomm = red->psubcomm->comm;
85: KSPCreate(subcomm,&red->ksp);
86: PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);
87: PetscLogObjectParent(pc,red->ksp);
88: KSPSetType(red->ksp,KSPPREONLY);
89: KSPGetPC(red->ksp,&red->pc);
90: PCSetType(red->pc,PCLU);
92: PCGetOptionsPrefix(pc,&prefix);
93: KSPSetOptionsPrefix(red->ksp,prefix);
94: KSPAppendOptionsPrefix(red->ksp,"redundant_");
95: } else {
96: subcomm = red->psubcomm->comm;
97: }
99: /* create working vectors xsub/ysub and xdup/ydup */
100: VecGetLocalSize(vec,&mlocal);
101: VecGetOwnershipRange(vec,&mstart,&mend);
103: /* get local size of xsub/ysub */
104: MPI_Comm_size(subcomm,&subsize);
105: MPI_Comm_rank(subcomm,&subrank);
106: rstart_sub = pc->pmat->rmap->range[red->psubcomm->n*subrank]; /* rstart in xsub/ysub */
107: if (subrank+1 < subsize){
108: rend_sub = pc->pmat->rmap->range[red->psubcomm->n*(subrank+1)];
109: } else {
110: rend_sub = m;
111: }
112: mloc_sub = rend_sub - rstart_sub;
113: VecCreateMPI(subcomm,mloc_sub,PETSC_DECIDE,&red->ysub);
114: /* create xsub with empty local arrays, because xdup's arrays will be placed into it */
115: VecCreateMPIWithArray(subcomm,mloc_sub,PETSC_DECIDE,PETSC_NULL,&red->xsub);
117: /* create xdup and ydup. ydup has empty local arrays because ysub's arrays will be place into it.
118: Note: we use communicator dupcomm, not ((PetscObject)pc)->comm! */
119: VecCreateMPI(red->psubcomm->dupparent,mloc_sub,PETSC_DECIDE,&red->xdup);
120: VecCreateMPIWithArray(red->psubcomm->dupparent,mloc_sub,PETSC_DECIDE,PETSC_NULL,&red->ydup);
121:
122: /* create vec scatters */
123: if (!red->scatterin){
124: IS is1,is2;
125: PetscInt *idx1,*idx2,i,j,k;
127: PetscMalloc(2*red->psubcomm->n*mlocal*sizeof(PetscInt),&idx1);
128: idx2 = idx1 + red->psubcomm->n*mlocal;
129: j = 0;
130: for (k=0; k<red->psubcomm->n; k++){
131: for (i=mstart; i<mend; i++){
132: idx1[j] = i;
133: idx2[j++] = i + m*k;
134: }
135: }
136: ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx1,&is1);
137: ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx2,&is2);
138: VecScatterCreate(vec,is1,red->xdup,is2,&red->scatterin);
139: ISDestroy(is1);
140: ISDestroy(is2);
142: ISCreateStride(comm,mlocal,mstart+ red->psubcomm->color*m,1,&is1);
143: ISCreateStride(comm,mlocal,mstart,1,&is2);
144: VecScatterCreate(red->xdup,is1,vec,is2,&red->scatterout);
145: ISDestroy(is1);
146: ISDestroy(is2);
147: PetscFree(idx1);
148: }
149: }
150: VecDestroy(vec);
152: /* if pmatrix set by user is sequential then we do not need to gather the parallel matrix */
153: MPI_Comm_size(comm,&size);
154: if (size == 1) {
155: red->useparallelmat = PETSC_FALSE;
156: }
158: if (red->useparallelmat) {
159: if (pc->setupcalled == 1 && pc->flag == DIFFERENT_NONZERO_PATTERN) {
160: /* destroy old matrices */
161: if (red->pmats) {
162: MatDestroy(red->pmats);
163: }
164: } else if (pc->setupcalled == 1) {
165: reuse = MAT_REUSE_MATRIX;
166: str = SAME_NONZERO_PATTERN;
167: }
168:
169: /* grab the parallel matrix and put it into processors of a subcomminicator */
170: /*--------------------------------------------------------------------------*/
171: VecGetLocalSize(red->ysub,&mlocal_sub);
172: MatGetRedundantMatrix(pc->pmat,red->psubcomm->n,red->psubcomm->comm,mlocal_sub,reuse,&red->pmats);
173: /* tell PC of the subcommunicator its operators */
174: KSPSetOperators(red->ksp,red->pmats,red->pmats,str);
175: } else {
176: KSPSetOperators(red->ksp,pc->mat,pc->pmat,pc->flag);
177: }
178: if (pc->setfromoptionscalled){
179: KSPSetFromOptions(red->ksp);
180: }
181: KSPSetUp(red->ksp);
182: return(0);
183: }
187: static PetscErrorCode PCApply_Redundant(PC pc,Vec x,Vec y)
188: {
189: PC_Redundant *red = (PC_Redundant*)pc->data;
191: PetscScalar *array;
194: /* scatter x to xdup */
195: VecScatterBegin(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);
196: VecScatterEnd(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);
197:
198: /* place xdup's local array into xsub */
199: VecGetArray(red->xdup,&array);
200: VecPlaceArray(red->xsub,(const PetscScalar*)array);
202: /* apply preconditioner on each processor */
203: PCApply(red->pc,red->xsub,red->ysub);
204: VecResetArray(red->xsub);
205: VecRestoreArray(red->xdup,&array);
206:
207: /* place ysub's local array into ydup */
208: VecGetArray(red->ysub,&array);
209: VecPlaceArray(red->ydup,(const PetscScalar*)array);
211: /* scatter ydup to y */
212: VecScatterBegin(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);
213: VecScatterEnd(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);
214: VecResetArray(red->ydup);
215: VecRestoreArray(red->ysub,&array);
216: return(0);
217: }
221: static PetscErrorCode PCDestroy_Redundant(PC pc)
222: {
223: PC_Redundant *red = (PC_Redundant*)pc->data;
227: if (red->scatterin) {VecScatterDestroy(red->scatterin);}
228: if (red->scatterout) {VecScatterDestroy(red->scatterout);}
229: if (red->ysub) {VecDestroy(red->ysub);}
230: if (red->xsub) {VecDestroy(red->xsub);}
231: if (red->xdup) {VecDestroy(red->xdup);}
232: if (red->ydup) {VecDestroy(red->ydup);}
233: if (red->pmats) {
234: MatDestroy(red->pmats);
235: }
236: if (red->psubcomm) {PetscSubcommDestroy(red->psubcomm);}
237: if (red->ksp) {KSPDestroy(red->ksp);}
238: PetscFree(red);
239: return(0);
240: }
244: static PetscErrorCode PCSetFromOptions_Redundant(PC pc)
245: {
247: PC_Redundant *red = (PC_Redundant*)pc->data;
250: PetscOptionsHead("Redundant options");
251: PetscOptionsInt("-pc_redundant_number","Number of redundant pc","PCRedundantSetNumber",red->nsubcomm,&red->nsubcomm,0);
252: PetscOptionsTail();
253: return(0);
254: }
259: PetscErrorCode PCRedundantSetNumber_Redundant(PC pc,PetscInt nreds)
260: {
261: PC_Redundant *red = (PC_Redundant*)pc->data;
264: red->nsubcomm = nreds;
265: return(0);
266: }
271: /*@
272: PCRedundantSetNumber - Sets the number of redundant preconditioner contexts.
274: Collective on PC
276: Input Parameters:
277: + pc - the preconditioner context
278: - nredundant - number of redundant preconditioner contexts; for example if you are using 64 MPI processes and
279: use an nredundant of 4 there will be 4 parallel solves each on 16 = 64/4 processes.
281: Level: advanced
283: .keywords: PC, redundant solve
284: @*/
285: PetscErrorCode PCRedundantSetNumber(PC pc,PetscInt nredundant)
286: {
287: PetscErrorCode ierr,(*f)(PC,PetscInt);
291: if (nredundant <= 0) SETERRQ1(PETSC_ERR_ARG_WRONG, "num of redundant pc %D must be positive",nredundant);
292: PetscObjectQueryFunction((PetscObject)pc,"PCRedundantSetNumber_C",(void (**)(void))&f);
293: if (f) {
294: (*f)(pc,nredundant);
295: }
296: return(0);
297: }
302: PetscErrorCode PCRedundantSetScatter_Redundant(PC pc,VecScatter in,VecScatter out)
303: {
304: PC_Redundant *red = (PC_Redundant*)pc->data;
308: PetscObjectReference((PetscObject)in);
309: if (red->scatterin) { VecScatterDestroy(red->scatterin); }
310: red->scatterin = in;
311: PetscObjectReference((PetscObject)out);
312: if (red->scatterout) { VecScatterDestroy(red->scatterout); }
313: red->scatterout = out;
314: return(0);
315: }
320: /*@
321: PCRedundantSetScatter - Sets the scatter used to copy values into the
322: redundant local solve and the scatter to move them back into the global
323: vector.
325: Collective on PC
327: Input Parameters:
328: + pc - the preconditioner context
329: . in - the scatter to move the values in
330: - out - the scatter to move them out
332: Level: advanced
334: .keywords: PC, redundant solve
335: @*/
336: PetscErrorCode PCRedundantSetScatter(PC pc,VecScatter in,VecScatter out)
337: {
338: PetscErrorCode ierr,(*f)(PC,VecScatter,VecScatter);
344: PetscObjectQueryFunction((PetscObject)pc,"PCRedundantSetScatter_C",(void (**)(void))&f);
345: if (f) {
346: (*f)(pc,in,out);
347: }
348: return(0);
349: }
354: PetscErrorCode PCRedundantGetPC_Redundant(PC pc,PC *innerpc)
355: {
357: PC_Redundant *red = (PC_Redundant*)pc->data;
358: MPI_Comm comm,subcomm;
359: const char *prefix;
362: if (!red->psubcomm) {
363: PetscObjectGetComm((PetscObject)pc,&comm);
364: PetscSubcommCreate(comm,red->nsubcomm,&red->psubcomm);
365: PetscLogObjectMemory(pc,sizeof(PetscSubcomm));
367: /* create a new PC that processors in each subcomm have copy of */
368: subcomm = red->psubcomm->comm;
369: KSPCreate(subcomm,&red->ksp);
370: PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);
371: PetscLogObjectParent(pc,red->ksp);
372: KSPSetType(red->ksp,KSPPREONLY);
373: KSPGetPC(red->ksp,&red->pc);
374: PCSetType(red->pc,PCLU);
376: PCGetOptionsPrefix(pc,&prefix);
377: KSPSetOptionsPrefix(red->ksp,prefix);
378: KSPAppendOptionsPrefix(red->ksp,"redundant_");
379: }
381: KSPGetPC(red->ksp,innerpc);
382: return(0);
383: }
388: /*@
389: PCRedundantGetPC - Gets the sequential PC created by the redundant PC.
391: Not Collective
393: Input Parameter:
394: . pc - the preconditioner context
396: Output Parameter:
397: . innerpc - the sequential PC
399: Level: advanced
401: .keywords: PC, redundant solve
402: @*/
403: PetscErrorCode PCRedundantGetPC(PC pc,PC *innerpc)
404: {
405: PetscErrorCode ierr,(*f)(PC,PC*);
410: PetscObjectQueryFunction((PetscObject)pc,"PCRedundantGetPC_C",(void (**)(void))&f);
411: if (f) {
412: (*f)(pc,innerpc);
413: }
414: return(0);
415: }
420: PetscErrorCode PCRedundantGetOperators_Redundant(PC pc,Mat *mat,Mat *pmat)
421: {
422: PC_Redundant *red = (PC_Redundant*)pc->data;
425: if (mat) *mat = red->pmats;
426: if (pmat) *pmat = red->pmats;
427: return(0);
428: }
433: /*@
434: PCRedundantGetOperators - gets the sequential matrix and preconditioner matrix
436: Not Collective
438: Input Parameter:
439: . pc - the preconditioner context
441: Output Parameters:
442: + mat - the matrix
443: - pmat - the (possibly different) preconditioner matrix
445: Level: advanced
447: .keywords: PC, redundant solve
448: @*/
449: PetscErrorCode PCRedundantGetOperators(PC pc,Mat *mat,Mat *pmat)
450: {
451: PetscErrorCode ierr,(*f)(PC,Mat*,Mat*);
457: PetscObjectQueryFunction((PetscObject)pc,"PCRedundantGetOperators_C",(void (**)(void))&f);
458: if (f) {
459: (*f)(pc,mat,pmat);
460: }
461: return(0);
462: }
464: /* -------------------------------------------------------------------------------------*/
465: /*MC
466: PCREDUNDANT - Runs a preconditioner for the entire problem on subgroups of processors
468: Options for the redundant preconditioners can be set with -redundant_pc_xxx
470: Options Database:
471: . -pc_redundant_number <n> - number of redundant solves, for example if you are using 64 MPI processes and
472: use an n of 4 there will be 4 parallel solves each on 16 = 64/4 processes.
474: Level: intermediate
476: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCRedundantSetScatter(),
477: PCRedundantGetPC(), PCRedundantGetOperators(), PCRedundantSetNumber()
478: M*/
483: PetscErrorCode PCCreate_Redundant(PC pc)
484: {
486: PC_Redundant *red;
487: PetscMPIInt size;
488:
490: PetscNewLog(pc,PC_Redundant,&red);
491: MPI_Comm_size(((PetscObject)pc)->comm,&size);
492: red->nsubcomm = size;
493: red->useparallelmat = PETSC_TRUE;
494: pc->data = (void*)red;
496: pc->ops->apply = PCApply_Redundant;
497: pc->ops->applytranspose = 0;
498: pc->ops->setup = PCSetUp_Redundant;
499: pc->ops->destroy = PCDestroy_Redundant;
500: pc->ops->setfromoptions = PCSetFromOptions_Redundant;
501: pc->ops->view = PCView_Redundant;
502: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantSetScatter_C","PCRedundantSetScatter_Redundant",
503: PCRedundantSetScatter_Redundant);
504: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantSetNumber_C","PCRedundantSetNumber_Redundant",
505: PCRedundantSetNumber_Redundant);
506: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantGetPC_C","PCRedundantGetPC_Redundant",
507: PCRedundantGetPC_Redundant);
508: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantGetOperators_C","PCRedundantGetOperators_Redundant",
509: PCRedundantGetOperators_Redundant);
510: return(0);
511: }