Actual source code: fieldsplit.c
1: /*
3: */
4: #include src/ksp/pc/pcimpl.h
6: typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
7: struct _PC_FieldSplitLink {
8: KSP ksp;
9: Vec x,y;
10: PetscInt nfields;
11: PetscInt *fields;
12: VecScatter sctx;
13: PC_FieldSplitLink next;
14: };
16: typedef struct {
17: PCCompositeType type; /* additive or multiplicative */
18: PetscTruth defaultsplit;
19: PetscInt bs;
20: PetscInt nsplits;
21: Vec *x,*y,w1,w2;
22: Mat *pmat;
23: IS *is;
24: PC_FieldSplitLink head;
25: } PC_FieldSplit;
29: static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
30: {
31: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
32: PetscErrorCode ierr;
33: PetscTruth iascii;
34: PetscInt i,j;
35: PC_FieldSplitLink link = jac->head;
36: const char *types[] = {"additive","multiplicative"};
39: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
40: if (iascii) {
41: PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D",types[jac->type],jac->nsplits);
42: PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");
43: PetscViewerASCIIPushTab(viewer);
44: for (i=0; i<jac->nsplits; i++) {
45: PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);
46: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
47: for (j=0; j<link->nfields; j++) {
48: if (j > 0) {
49: PetscViewerASCIIPrintf(viewer,",");
50: }
51: PetscViewerASCIIPrintf(viewer," %D",link->fields[j]);
52: }
53: PetscViewerASCIIPrintf(viewer,"\n");
54: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
55: KSPView(link->ksp,viewer);
56: link = link->next;
57: }
58: PetscViewerASCIIPopTab(viewer);
59: } else {
60: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
61: }
62: return(0);
63: }
67: static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
68: {
69: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
70: PetscErrorCode ierr;
71: PC_FieldSplitLink link = jac->head;
72: PetscInt i;
73: PetscTruth flg = PETSC_FALSE,*fields;
76: PetscOptionsGetLogical(pc->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);
77: if (!link || flg) {
78: PetscLogInfo(pc,"PCFieldSplitSetDefaults: Using default splitting of fields");
79: if (jac->bs <= 0) {
80: MatGetBlockSize(pc->pmat,&jac->bs);
81: }
82: PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);
83: PetscMemzero(fields,jac->bs*sizeof(PetscTruth));
84: while (link) {
85: for (i=0; i<link->nfields; i++) {
86: fields[link->fields[i]] = PETSC_TRUE;
87: }
88: link = link->next;
89: }
90: jac->defaultsplit = PETSC_TRUE;
91: for (i=0; i<jac->bs; i++) {
92: if (!fields[i]) {
93: PCFieldSplitSetFields(pc,1,&i);
94: } else {
95: jac->defaultsplit = PETSC_FALSE;
96: }
97: }
98: }
99: return(0);
100: }
105: static PetscErrorCode PCSetUp_FieldSplit(PC pc)
106: {
107: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
108: PetscErrorCode ierr;
109: PC_FieldSplitLink link;
110: PetscInt i,nsplit;
111: MatStructure flag = pc->flag;
114: PCFieldSplitSetDefaults(pc);
115: nsplit = jac->nsplits;
116: link = jac->head;
118: /* get the matrices for each split */
119: if (!jac->is) {
120: PetscInt rstart,rend,nslots,bs;
122: MatGetBlockSize(pc->pmat,&bs);
123: MatGetOwnershipRange(pc->pmat,&rstart,&rend);
124: nslots = (rend - rstart)/bs;
125: PetscMalloc(nsplit*sizeof(IS),&jac->is);
126: for (i=0; i<nsplit; i++) {
127: if (jac->defaultsplit) {
128: ISCreateStride(pc->comm,nslots,rstart+i,nsplit,&jac->is[i]);
129: } else {
130: PetscInt *ii,j,k,nfields = link->nfields,*fields = link->fields;
131: PetscTruth sorted;
132: PetscMalloc(link->nfields*nslots*sizeof(PetscInt),&ii);
133: for (j=0; j<nslots; j++) {
134: for (k=0; k<nfields; k++) {
135: ii[nfields*j + k] = rstart + bs*j + fields[k];
136: }
137: }
138: ISCreateGeneral(pc->comm,nslots*nfields,ii,&jac->is[i]);
139: ISSorted(jac->is[i],&sorted);
140: if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split");
141: PetscFree(ii);
142: link = link->next;
143: }
144: }
145: }
146:
147: if (!jac->pmat) {
148: MatGetSubMatrices(pc->pmat,nsplit,jac->is,jac->is,MAT_INITIAL_MATRIX,&jac->pmat);
149: } else {
150: MatGetSubMatrices(pc->pmat,nsplit,jac->is,jac->is,MAT_REUSE_MATRIX,&jac->pmat);
151: }
153: /* set up the individual PCs */
154: i = 0;
155: link = jac->head;
156: while (link) {
157: KSPSetOperators(link->ksp,jac->pmat[i],jac->pmat[i],flag);
158: KSPSetFromOptions(link->ksp);
159: KSPSetUp(link->ksp);
160: i++;
161: link = link->next;
162: }
163:
164: /* create work vectors for each split */
165: if (!jac->x) {
166: Vec xtmp;
167: PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);
168: link = jac->head;
169: for (i=0; i<nsplit; i++) {
170: Mat A;
171: KSPGetOperators(link->ksp,PETSC_NULL,&A,PETSC_NULL);
172: MatGetVecs(A,&link->x,&link->y);
173: jac->x[i] = link->x;
174: jac->y[i] = link->y;
175: link = link->next;
176: }
177: /* compute scatter contexts needed by multiplicative versions and non-default splits */
178:
179: link = jac->head;
180: MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);
181: for (i=0; i<nsplit; i++) {
182: VecScatterCreate(xtmp,jac->is[i],jac->x[i],PETSC_NULL,&link->sctx);
183: link = link->next;
184: }
185: VecDestroy(xtmp);
186: }
187: return(0);
188: }
190: #define FieldSplitSplitSolveAdd(link,xx,yy) \
191: (VecScatterBegin(xx,link->x,INSERT_VALUES,SCATTER_FORWARD,link->sctx) || \
192: VecScatterEnd(xx,link->x,INSERT_VALUES,SCATTER_FORWARD,link->sctx) || \
193: KSPSolve(link->ksp,link->x,link->y) || \
194: VecScatterBegin(link->y,yy,ADD_VALUES,SCATTER_REVERSE,link->sctx) || \
195: VecScatterEnd(link->y,yy,ADD_VALUES,SCATTER_REVERSE,link->sctx))
199: static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
200: {
201: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
202: PetscErrorCode ierr;
203: PC_FieldSplitLink link = jac->head;
204: PetscScalar zero = 0.0,mone = -1.0;
207: if (jac->type == PC_COMPOSITE_ADDITIVE) {
208: if (jac->defaultsplit) {
209: VecStrideGatherAll(x,jac->x,INSERT_VALUES);
210: while (link) {
211: KSPSolve(link->ksp,link->x,link->y);
212: link = link->next;
213: }
214: VecStrideScatterAll(jac->y,y,INSERT_VALUES);
215: } else {
216: PetscInt i = 0;
218: VecSet(&zero,y);
219: while (link) {
220: FieldSplitSplitSolveAdd(link,x,y);
221: link = link->next;
222: i++;
223: }
224: }
225: } else {
226: if (!jac->w1) {
227: VecDuplicate(x,&jac->w1);
228: VecDuplicate(x,&jac->w2);
229: }
230: VecSet(&zero,y);
231: FieldSplitSplitSolveAdd(link,x,y);
232: while (link->next) {
233: link = link->next;
234: MatMult(pc->pmat,y,jac->w1);
235: VecWAXPY(&mone,jac->w1,x,jac->w2);
236: FieldSplitSplitSolveAdd(link,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 link = jac->head,next;
251: while (link) {
252: KSPDestroy(link->ksp);
253: if (link->x) {VecDestroy(link->x);}
254: if (link->y) {VecDestroy(link->y);}
255: if (link->sctx) {VecScatterDestroy(link->sctx);}
256: next = link->next;
257: PetscFree2(link,link->fields);
258: link = next;
259: }
260: if (jac->x) {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],indx;
280: PetscTruth flg;
281: char optionname[128];
282: const char *types[] = {"additive","multiplicative"};
285: PetscOptionsHead("FieldSplit options");
286: PetscOptionsEList("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",types,2,types[0],&indx,&flg);
287: if (flg) {
288: PCFieldSplitSetType(pc,(PCCompositeType)indx);
289: }
290: while (PETSC_TRUE) {
291: sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
292: nfields = 12;
293: PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);
294: if (!flg) break;
295: if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields");
296: PCFieldSplitSetFields(pc,nfields,fields);
297: }
298: PetscOptionsTail();
299: return(0);
300: }
302: /*------------------------------------------------------------------------------------*/
307: PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields)
308: {
309: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
310: PetscErrorCode ierr;
311: PC_FieldSplitLink link,next = jac->head;
312: char prefix[128];
315: if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested");
316: PetscMalloc2(1,struct _PC_FieldSplitLink,&link,n,PetscInt,&link->fields);
317: PetscMemcpy(link->fields,fields,n*sizeof(PetscInt));
318: link->nfields = n;
319: link->next = PETSC_NULL;
320: KSPCreate(pc->comm,&link->ksp);
321: KSPSetType(link->ksp,KSPPREONLY);
323: if (pc->prefix) {
324: sprintf(prefix,"%sfieldsplit_%d_",pc->prefix,(int)jac->nsplits);
325: } else {
326: sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
327: }
328: KSPSetOptionsPrefix(link->ksp,prefix);
330: if (!next) {
331: jac->head = link;
332: } else {
333: while (next->next) {
334: next = next->next;
335: }
336: next->next = link;
337: }
338: jac->nsplits++;
339: return(0);
340: }
347: PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
348: {
349: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
350: PetscErrorCode ierr;
351: PetscInt cnt = 0;
352: PC_FieldSplitLink link = jac->head;
355: PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);
356: while (link) {
357: (*subksp)[cnt++] = link->ksp;
358: link = link->next;
359: }
360: if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits);
361: *n = jac->nsplits;
362: return(0);
363: }
368: /*@
369: PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
371: Collective on PC
373: Input Parameters:
374: + pc - the preconditioner context
375: . n - the number of fields in this split
376: . fields - the fields in this split
378: Level: intermediate
380: .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT
382: @*/
383: PetscErrorCode PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields)
384: {
385: PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *);
389: PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);
390: if (f) {
391: (*f)(pc,n,fields);
392: }
393: return(0);
394: }
398: /*@C
399: PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
400:
401: Collective on KSP
403: Input Parameter:
404: . pc - the preconditioner context
406: Output Parameters:
407: + n - the number of split
408: - pc - the array of KSP contexts
410: Note:
411: After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed
413: You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
415: Level: advanced
417: .seealso: PCFIELDSPLIT
418: @*/
419: PetscErrorCode PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
420: {
421: PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **);
426: PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);
427: if (f) {
428: (*f)(pc,n,subksp);
429: } else {
430: SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC");
431: }
432: return(0);
433: }
438: PetscErrorCode PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
439: {
440: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
443: jac->type = type;
444: return(0);
445: }
450: /*@C
451: PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
452:
453: Collective on PC
455: Input Parameter:
456: . pc - the preconditioner context
457: . type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE
459: Options Database Key:
460: . -pc_fieldsplit_type <type: one of multiplicative, additive, special> - Sets fieldsplit preconditioner type
462: Level: Developer
464: .keywords: PC, set, type, composite preconditioner, additive, multiplicative
466: .seealso: PCCompositeSetType()
468: @*/
469: PetscErrorCode PCFieldSplitSetType(PC pc,PCCompositeType type)
470: {
471: PetscErrorCode ierr,(*f)(PC,PCCompositeType);
475: PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);
476: if (f) {
477: (*f)(pc,type);
478: }
479: return(0);
480: }
482: /* -------------------------------------------------------------------------------------*/
483: /*MC
484: PCFieldSplit - Preconditioner created by combining seperate preconditioners for individual
485: fields or groups of fields
488: To set options on the solvers for each block append -sub_ to all the PC
489: options database keys. For example, -sub_pc_type ilu -sub_pc_ilu_levels 1
490:
491: To set the options on the solvers seperate for each block call PCFieldSplitGetSubKSP()
492: and set the options directly on the resulting KSP object
494: Level: intermediate
496: Options Database Keys:
497: + -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
498: . -pc_splitfield_default - automatically add any fields to additional splits that have not
499: been supplied explicitly by -pc_splitfield_%d_fields
500: - -pc_splitfield_type <additive,multiplicative>
502: Concepts: physics based preconditioners
504: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC,
505: PCFieldSplitGetSubKSP(), PCFieldSplitSetFields()
506: M*/
511: PetscErrorCode PCCreate_FieldSplit(PC pc)
512: {
514: PC_FieldSplit *jac;
517: PetscNew(PC_FieldSplit,&jac);
518: PetscLogObjectMemory(pc,sizeof(PC_FieldSplit));
519: PetscMemzero(jac,sizeof(PC_FieldSplit));
520: jac->bs = -1;
521: jac->nsplits = 0;
522: jac->type = PC_COMPOSITE_ADDITIVE;
523: pc->data = (void*)jac;
525: pc->ops->apply = PCApply_FieldSplit;
526: pc->ops->setup = PCSetUp_FieldSplit;
527: pc->ops->destroy = PCDestroy_FieldSplit;
528: pc->ops->setfromoptions = PCSetFromOptions_FieldSplit;
529: pc->ops->view = PCView_FieldSplit;
530: pc->ops->applyrichardson = 0;
532: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
533: PCFieldSplitGetSubKSP_FieldSplit);
534: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
535: PCFieldSplitSetFields_FieldSplit);
536: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
537: PCFieldSplitSetType_FieldSplit);
538: return(0);
539: }