Actual source code: fieldsplit.c
1: #define PETSCKSP_DLL
3: /*
5: */
6: #include private/pcimpl.h
8: typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
9: struct _PC_FieldSplitLink {
10: KSP ksp;
11: Vec x,y;
12: PetscInt nfields;
13: PetscInt *fields;
14: VecScatter sctx;
15: PC_FieldSplitLink next;
16: };
18: typedef struct {
19: PCCompositeType type; /* additive or multiplicative */
20: PetscTruth defaultsplit;
21: PetscInt bs;
22: PetscInt nsplits;
23: Vec *x,*y,w1,w2;
24: Mat *pmat;
25: IS *is;
26: PC_FieldSplitLink head;
27: } PC_FieldSplit;
31: static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
32: {
33: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
34: PetscErrorCode ierr;
35: PetscTruth iascii;
36: PetscInt i,j;
37: PC_FieldSplitLink ilink = jac->head;
40: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
41: if (iascii) {
42: PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D",PCCompositeTypes[jac->type],jac->nsplits);
43: PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");
44: PetscViewerASCIIPushTab(viewer);
45: for (i=0; i<jac->nsplits; i++) {
46: PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);
47: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
48: for (j=0; j<ilink->nfields; j++) {
49: if (j > 0) {
50: PetscViewerASCIIPrintf(viewer,",");
51: }
52: PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);
53: }
54: PetscViewerASCIIPrintf(viewer,"\n");
55: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
56: KSPView(ilink->ksp,viewer);
57: ilink = ilink->next;
58: }
59: PetscViewerASCIIPopTab(viewer);
60: } else {
61: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
62: }
63: return(0);
64: }
68: static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
69: {
70: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
71: PetscErrorCode ierr;
72: PC_FieldSplitLink ilink = jac->head;
73: PetscInt i;
74: PetscTruth flg = PETSC_FALSE,*fields;
77: PetscOptionsGetTruth(pc->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);
78: if (!ilink || flg) {
79: PetscInfo(pc,"Using default splitting of fields\n");
80: if (jac->bs <= 0) {
81: MatGetBlockSize(pc->pmat,&jac->bs);
82: }
83: PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);
84: PetscMemzero(fields,jac->bs*sizeof(PetscTruth));
85: while (ilink) {
86: for (i=0; i<ilink->nfields; i++) {
87: fields[ilink->fields[i]] = PETSC_TRUE;
88: }
89: ilink = ilink->next;
90: }
91: jac->defaultsplit = PETSC_TRUE;
92: for (i=0; i<jac->bs; i++) {
93: if (!fields[i]) {
94: PCFieldSplitSetFields(pc,1,&i);
95: } else {
96: jac->defaultsplit = PETSC_FALSE;
97: }
98: }
99: }
100: return(0);
101: }
106: static PetscErrorCode PCSetUp_FieldSplit(PC pc)
107: {
108: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
109: PetscErrorCode ierr;
110: PC_FieldSplitLink ilink;
111: PetscInt i,nsplit;
112: MatStructure flag = pc->flag;
115: PCFieldSplitSetDefaults(pc);
116: nsplit = jac->nsplits;
117: ilink = jac->head;
119: /* get the matrices for each split */
120: if (!jac->is) {
121: PetscInt rstart,rend,nslots,bs;
123: MatGetBlockSize(pc->pmat,&bs);
124: MatGetOwnershipRange(pc->pmat,&rstart,&rend);
125: nslots = (rend - rstart)/bs;
126: PetscMalloc(nsplit*sizeof(IS),&jac->is);
127: for (i=0; i<nsplit; i++) {
128: if (jac->defaultsplit) {
129: ISCreateStride(pc->comm,nslots,rstart+i,nsplit,&jac->is[i]);
130: } else {
131: PetscInt *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
132: PetscTruth sorted;
133: PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);
134: for (j=0; j<nslots; j++) {
135: for (k=0; k<nfields; k++) {
136: ii[nfields*j + k] = rstart + bs*j + fields[k];
137: }
138: }
139: ISCreateGeneral(pc->comm,nslots*nfields,ii,&jac->is[i]);
140: ISSorted(jac->is[i],&sorted);
141: if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split");
142: PetscFree(ii);
143: ilink = ilink->next;
144: }
145: }
146: }
147:
148: if (!jac->pmat) {
149: MatGetSubMatrices(pc->pmat,nsplit,jac->is,jac->is,MAT_INITIAL_MATRIX,&jac->pmat);
150: } else {
151: MatGetSubMatrices(pc->pmat,nsplit,jac->is,jac->is,MAT_REUSE_MATRIX,&jac->pmat);
152: }
154: /* set up the individual PCs */
155: i = 0;
156: ilink = jac->head;
157: while (ilink) {
158: KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);
159: KSPSetFromOptions(ilink->ksp);
160: KSPSetUp(ilink->ksp);
161: i++;
162: ilink = ilink->next;
163: }
164:
165: /* create work vectors for each split */
166: if (!jac->x) {
167: Vec xtmp;
168: PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);
169: ilink = jac->head;
170: for (i=0; i<nsplit; i++) {
171: Mat A;
172: KSPGetOperators(ilink->ksp,PETSC_NULL,&A,PETSC_NULL);
173: MatGetVecs(A,&ilink->x,&ilink->y);
174: jac->x[i] = ilink->x;
175: jac->y[i] = ilink->y;
176: ilink = ilink->next;
177: }
178: /* compute scatter contexts needed by multiplicative versions and non-default splits */
179:
180: ilink = jac->head;
181: MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);
182: for (i=0; i<nsplit; i++) {
183: VecScatterCreate(xtmp,jac->is[i],jac->x[i],PETSC_NULL,&ilink->sctx);
184: ilink = ilink->next;
185: }
186: VecDestroy(xtmp);
187: }
188: return(0);
189: }
191: #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
192: (VecScatterBegin(xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \
193: VecScatterEnd(xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \
194: KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
195: VecScatterBegin(ilink->y,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx) || \
196: VecScatterEnd(ilink->y,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx))
200: static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
201: {
202: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
203: PetscErrorCode ierr;
204: PC_FieldSplitLink ilink = jac->head;
207: if (jac->type == PC_COMPOSITE_ADDITIVE) {
208: if (jac->defaultsplit) {
209: VecStrideGatherAll(x,jac->x,INSERT_VALUES);
210: while (ilink) {
211: KSPSolve(ilink->ksp,ilink->x,ilink->y);
212: ilink = ilink->next;
213: }
214: VecStrideScatterAll(jac->y,y,INSERT_VALUES);
215: } else {
216: PetscInt i = 0;
218: VecSet(y,0.0);
219: while (ilink) {
220: FieldSplitSplitSolveAdd(ilink,x,y);
221: ilink = ilink->next;
222: i++;
223: }
224: }
225: } else {
226: if (!jac->w1) {
227: VecDuplicate(x,&jac->w1);
228: VecDuplicate(x,&jac->w2);
229: }
230: VecSet(y,0.0);
231: FieldSplitSplitSolveAdd(ilink,x,y);
232: while (ilink->next) {
233: ilink = ilink->next;
234: MatMult(pc->pmat,y,jac->w1);
235: VecWAXPY(jac->w2,-1.0,jac->w1,x);
236: FieldSplitSplitSolveAdd(ilink,jac->w2,y);
237: }
238: }
239: return(0);
240: }
244: static PetscErrorCode PCDestroy_FieldSplit(PC pc)
245: {
246: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
247: PetscErrorCode ierr;
248: PC_FieldSplitLink ilink = jac->head,next;
251: while (ilink) {
252: KSPDestroy(ilink->ksp);
253: if (ilink->x) {VecDestroy(ilink->x);}
254: if (ilink->y) {VecDestroy(ilink->y);}
255: if (ilink->sctx) {VecScatterDestroy(ilink->sctx);}
256: next = ilink->next;
257: PetscFree2(ilink,ilink->fields);
258: ilink = next;
259: }
260: PetscFree2(jac->x,jac->y);
261: if (jac->pmat) {MatDestroyMatrices(jac->nsplits,&jac->pmat);}
262: if (jac->is) {
263: PetscInt i;
264: for (i=0; i<jac->nsplits; i++) {ISDestroy(jac->is[i]);}
265: PetscFree(jac->is);
266: }
267: if (jac->w1) {VecDestroy(jac->w1);}
268: if (jac->w2) {VecDestroy(jac->w2);}
269: PetscFree(jac);
270: return(0);
271: }
275: static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
276: /* This does not call KSPSetFromOptions() on the subksp's, see PCSetFromOptionsBJacobi/ASM() */
277: {
279: PetscInt i = 0,nfields,fields[12];
280: PetscTruth flg;
281: char optionname[128];
282: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
285: PetscOptionsHead("FieldSplit options");
286: PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);
287: while (PETSC_TRUE) {
288: sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
289: nfields = 12;
290: PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);
291: if (!flg) break;
292: if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields");
293: PCFieldSplitSetFields(pc,nfields,fields);
294: }
295: PetscOptionsTail();
296: return(0);
297: }
299: /*------------------------------------------------------------------------------------*/
304: PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields)
305: {
306: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
307: PetscErrorCode ierr;
308: PC_FieldSplitLink ilink,next = jac->head;
309: char prefix[128];
312: if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested");
313: PetscMalloc2(1,struct _PC_FieldSplitLink,&ilink,n,PetscInt,&ilink->fields);
314: PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));
315: ilink->nfields = n;
316: ilink->next = PETSC_NULL;
317: KSPCreate(pc->comm,&ilink->ksp);
318: KSPSetType(ilink->ksp,KSPPREONLY);
320: if (pc->prefix) {
321: sprintf(prefix,"%sfieldsplit_%d_",pc->prefix,(int)jac->nsplits);
322: } else {
323: sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
324: }
325: KSPSetOptionsPrefix(ilink->ksp,prefix);
327: if (!next) {
328: jac->head = ilink;
329: } else {
330: while (next->next) {
331: next = next->next;
332: }
333: next->next = ilink;
334: }
335: jac->nsplits++;
336: return(0);
337: }
344: PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
345: {
346: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
347: PetscErrorCode ierr;
348: PetscInt cnt = 0;
349: PC_FieldSplitLink ilink = jac->head;
352: PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);
353: while (ilink) {
354: (*subksp)[cnt++] = ilink->ksp;
355: ilink = ilink->next;
356: }
357: if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits);
358: *n = jac->nsplits;
359: return(0);
360: }
365: /*@
366: PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
368: Collective on PC
370: Input Parameters:
371: + pc - the preconditioner context
372: . n - the number of fields in this split
373: . fields - the fields in this split
375: Level: intermediate
377: .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT
379: @*/
380: PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields)
381: {
382: PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *);
386: PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);
387: if (f) {
388: (*f)(pc,n,fields);
389: }
390: return(0);
391: }
395: /*@C
396: PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
397:
398: Collective on KSP
400: Input Parameter:
401: . pc - the preconditioner context
403: Output Parameters:
404: + n - the number of split
405: - pc - the array of KSP contexts
407: Note:
408: After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed
410: You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
412: Level: advanced
414: .seealso: PCFIELDSPLIT
415: @*/
416: PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
417: {
418: PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **);
423: PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);
424: if (f) {
425: (*f)(pc,n,subksp);
426: } else {
427: SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC");
428: }
429: return(0);
430: }
435: PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
436: {
437: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
440: jac->type = type;
441: return(0);
442: }
447: /*@C
448: PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
449:
450: Collective on PC
452: Input Parameter:
453: . pc - the preconditioner context
454: . type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE
456: Options Database Key:
457: . -pc_fieldsplit_type <type: one of multiplicative, additive, special> - Sets fieldsplit preconditioner type
459: Level: Developer
461: .keywords: PC, set, type, composite preconditioner, additive, multiplicative
463: .seealso: PCCompositeSetType()
465: @*/
466: PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type)
467: {
468: PetscErrorCode ierr,(*f)(PC,PCCompositeType);
472: PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);
473: if (f) {
474: (*f)(pc,type);
475: }
476: return(0);
477: }
479: /* -------------------------------------------------------------------------------------*/
480: /*MC
481: PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
482: fields or groups of fields
485: To set options on the solvers for each block append -sub_ to all the PC
486: options database keys. For example, -sub_pc_type ilu -sub_pc_ilu_levels 1
487:
488: To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
489: and set the options directly on the resulting KSP object
491: Level: intermediate
493: Options Database Keys:
494: + -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
495: . -pc_splitfield_default - automatically add any fields to additional splits that have not
496: been supplied explicitly by -pc_splitfield_%d_fields
497: - -pc_splitfield_type <additive,multiplicative>
499: Concepts: physics based preconditioners
501: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC,
502: PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(),PCFieldSplitSetType()
503: M*/
508: PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc)
509: {
511: PC_FieldSplit *jac;
514: PetscNew(PC_FieldSplit,&jac);
515: PetscLogObjectMemory(pc,sizeof(PC_FieldSplit));
516: jac->bs = -1;
517: jac->nsplits = 0;
518: jac->type = PC_COMPOSITE_ADDITIVE;
519: pc->data = (void*)jac;
521: pc->ops->apply = PCApply_FieldSplit;
522: pc->ops->setup = PCSetUp_FieldSplit;
523: pc->ops->destroy = PCDestroy_FieldSplit;
524: pc->ops->setfromoptions = PCSetFromOptions_FieldSplit;
525: pc->ops->view = PCView_FieldSplit;
526: pc->ops->applyrichardson = 0;
528: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
529: PCFieldSplitGetSubKSP_FieldSplit);
530: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
531: PCFieldSplitSetFields_FieldSplit);
532: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
533: PCFieldSplitSetType_FieldSplit);
534: return(0);
535: }