Actual source code: fieldsplit.c
petsc-3.12.5 2020-03-29
2: #include <petsc/private/pcimpl.h>
3: #include <petsc/private/kspimpl.h> /* This is needed to provide the appropriate PETSC_EXTERN for KSP_Solve_FS ....*/
4: #include <petscdm.h>
6: const char *const PCFieldSplitSchurPreTypes[] = {"SELF","SELFP","A11","USER","FULL","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0};
7: const char *const PCFieldSplitSchurFactTypes[] = {"DIAG","LOWER","UPPER","FULL","PCFieldSplitSchurFactType","PC_FIELDSPLIT_SCHUR_FACT_",0};
9: PetscLogEvent KSP_Solve_FS_0,KSP_Solve_FS_1,KSP_Solve_FS_S,KSP_Solve_FS_U,KSP_Solve_FS_L,KSP_Solve_FS_2,KSP_Solve_FS_3,KSP_Solve_FS_4;
11: typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
12: struct _PC_FieldSplitLink {
13: KSP ksp;
14: Vec x,y,z;
15: char *splitname;
16: PetscInt nfields;
17: PetscInt *fields,*fields_col;
18: VecScatter sctx;
19: IS is,is_col;
20: PC_FieldSplitLink next,previous;
21: PetscLogEvent event;
22: };
24: typedef struct {
25: PCCompositeType type;
26: PetscBool defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */
27: PetscBool splitdefined; /* Flag is set after the splits have been defined, to prevent more splits from being added */
28: PetscInt bs; /* Block size for IS and Mat structures */
29: PetscInt nsplits; /* Number of field divisions defined */
30: Vec *x,*y,w1,w2;
31: Mat *mat; /* The diagonal block for each split */
32: Mat *pmat; /* The preconditioning diagonal block for each split */
33: Mat *Afield; /* The rows of the matrix associated with each split */
34: PetscBool issetup;
36: /* Only used when Schur complement preconditioning is used */
37: Mat B; /* The (0,1) block */
38: Mat C; /* The (1,0) block */
39: Mat schur; /* The Schur complement S = A11 - A10 A00^{-1} A01, the KSP here, kspinner, is H_1 in [El08] */
40: Mat schurp; /* Assembled approximation to S built by MatSchurComplement to be used as a preconditioning matrix when solving with S */
41: Mat schur_user; /* User-provided preconditioning matrix for the Schur complement */
42: PCFieldSplitSchurPreType schurpre; /* Determines which preconditioning matrix is used for the Schur complement */
43: PCFieldSplitSchurFactType schurfactorization;
44: KSP kspschur; /* The solver for S */
45: KSP kspupper; /* The solver for A in the upper diagonal part of the factorization (H_2 in [El08]) */
46: PetscScalar schurscale; /* Scaling factor for the Schur complement solution with DIAG factorization */
48: /* Only used when Golub-Kahan bidiagonalization preconditioning is used */
49: Mat H; /* The modified matrix H = A00 + nu*A01*A01' */
50: PetscReal gkbtol; /* Stopping tolerance for lower bound estimate */
51: PetscInt gkbdelay; /* The delay window for the stopping criterion */
52: PetscReal gkbnu; /* Parameter for augmented Lagrangian H = A + nu*A01*A01' */
53: PetscInt gkbmaxit; /* Maximum number of iterations for outer loop */
54: PetscBool gkbmonitor; /* Monitor for gkb iterations and the lower bound error */
55: PetscViewer gkbviewer; /* Viewer context for gkbmonitor */
56: Vec u,v,d,Hu; /* Work vectors for the GKB algorithm */
57: PetscScalar *vecz; /* Contains intermediate values, eg for lower bound */
59: PC_FieldSplitLink head;
60: PetscBool isrestrict; /* indicates PCFieldSplitRestrictIS() has been last called on this object, hack */
61: PetscBool suboptionsset; /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */
62: PetscBool dm_splits; /* Whether to use DMCreateFieldDecomposition() whenever possible */
63: PetscBool diag_use_amat; /* Whether to extract diagonal matrix blocks from Amat, rather than Pmat (weaker than -pc_use_amat) */
64: PetscBool offdiag_use_amat; /* Whether to extract off-diagonal matrix blocks from Amat, rather than Pmat (weaker than -pc_use_amat) */
65: PetscBool detect; /* Whether to form 2-way split by finding zero diagonal entries */
66: } PC_FieldSplit;
68: /*
69: Notes:
70: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
71: inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
72: PC you could change this.
73: */
75: /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it. This way the
76: * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */
77: static Mat FieldSplitSchurPre(PC_FieldSplit *jac)
78: {
79: switch (jac->schurpre) {
80: case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur;
81: case PC_FIELDSPLIT_SCHUR_PRE_SELFP: return jac->schurp;
82: case PC_FIELDSPLIT_SCHUR_PRE_A11: return jac->pmat[1];
83: case PC_FIELDSPLIT_SCHUR_PRE_FULL: /* We calculate this and store it in schur_user */
84: case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */
85: default:
86: return jac->schur_user ? jac->schur_user : jac->pmat[1];
87: }
88: }
91: #include <petscdraw.h>
92: static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
93: {
94: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
95: PetscErrorCode ierr;
96: PetscBool iascii,isdraw;
97: PetscInt i,j;
98: PC_FieldSplitLink ilink = jac->head;
101: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
102: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
103: if (iascii) {
104: if (jac->bs > 0) {
105: PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);
106: } else {
107: PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D\n",PCCompositeTypes[jac->type],jac->nsplits);
108: }
109: if (pc->useAmat) {
110: PetscViewerASCIIPrintf(viewer," using Amat (not Pmat) as operator for blocks\n");
111: }
112: if (jac->diag_use_amat) {
113: PetscViewerASCIIPrintf(viewer," using Amat (not Pmat) as operator for diagonal blocks\n");
114: }
115: if (jac->offdiag_use_amat) {
116: PetscViewerASCIIPrintf(viewer," using Amat (not Pmat) as operator for off-diagonal blocks\n");
117: }
118: PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");
119: for (i=0; i<jac->nsplits; i++) {
120: if (ilink->fields) {
121: PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);
122: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
123: for (j=0; j<ilink->nfields; j++) {
124: if (j > 0) {
125: PetscViewerASCIIPrintf(viewer,",");
126: }
127: PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);
128: }
129: PetscViewerASCIIPrintf(viewer,"\n");
130: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
131: } else {
132: PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);
133: }
134: KSPView(ilink->ksp,viewer);
135: ilink = ilink->next;
136: }
137: }
139: if (isdraw) {
140: PetscDraw draw;
141: PetscReal x,y,w,wd;
143: PetscViewerDrawGetDraw(viewer,0,&draw);
144: PetscDrawGetCurrentPoint(draw,&x,&y);
145: w = 2*PetscMin(1.0 - x,x);
146: wd = w/(jac->nsplits + 1);
147: x = x - wd*(jac->nsplits-1)/2.0;
148: for (i=0; i<jac->nsplits; i++) {
149: PetscDrawPushCurrentPoint(draw,x,y);
150: KSPView(ilink->ksp,viewer);
151: PetscDrawPopCurrentPoint(draw);
152: x += wd;
153: ilink = ilink->next;
154: }
155: }
156: return(0);
157: }
159: static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
160: {
161: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
162: PetscErrorCode ierr;
163: PetscBool iascii,isdraw;
164: PetscInt i,j;
165: PC_FieldSplitLink ilink = jac->head;
166: MatSchurComplementAinvType atype;
169: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
170: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
171: if (iascii) {
172: if (jac->bs > 0) {
173: PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactTypes[jac->schurfactorization]);
174: } else {
175: PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, factorization %s\n",PCFieldSplitSchurFactTypes[jac->schurfactorization]);
176: }
177: if (pc->useAmat) {
178: PetscViewerASCIIPrintf(viewer," using Amat (not Pmat) as operator for blocks\n");
179: }
180: switch (jac->schurpre) {
181: case PC_FIELDSPLIT_SCHUR_PRE_SELF:
182: PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from S itself\n");
183: break;
184: case PC_FIELDSPLIT_SCHUR_PRE_SELFP:
185: MatSchurComplementGetAinvType(jac->schur,&atype);
186: PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from Sp, an assembled approximation to S, which uses A00's %sdiagonal's inverse\n",atype == MAT_SCHUR_COMPLEMENT_AINV_DIAG ? "" : (atype == MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG ? "block " : "lumped "));break;
187: case PC_FIELDSPLIT_SCHUR_PRE_A11:
188: PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from A11\n");
189: break;
190: case PC_FIELDSPLIT_SCHUR_PRE_FULL:
191: PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from the exact Schur complement\n");
192: break;
193: case PC_FIELDSPLIT_SCHUR_PRE_USER:
194: if (jac->schur_user) {
195: PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from user provided matrix\n");
196: } else {
197: PetscViewerASCIIPrintf(viewer," Preconditioner for the Schur complement formed from A11\n");
198: }
199: break;
200: default:
201: SETERRQ1(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Invalid Schur preconditioning type: %d", jac->schurpre);
202: }
203: PetscViewerASCIIPrintf(viewer," Split info:\n");
204: PetscViewerASCIIPushTab(viewer);
205: for (i=0; i<jac->nsplits; i++) {
206: if (ilink->fields) {
207: PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);
208: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
209: for (j=0; j<ilink->nfields; j++) {
210: if (j > 0) {
211: PetscViewerASCIIPrintf(viewer,",");
212: }
213: PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);
214: }
215: PetscViewerASCIIPrintf(viewer,"\n");
216: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
217: } else {
218: PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);
219: }
220: ilink = ilink->next;
221: }
222: PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block\n");
223: PetscViewerASCIIPushTab(viewer);
224: if (jac->head) {
225: KSPView(jac->head->ksp,viewer);
226: } else {PetscViewerASCIIPrintf(viewer," not yet available\n");}
227: PetscViewerASCIIPopTab(viewer);
228: if (jac->head && jac->kspupper != jac->head->ksp) {
229: PetscViewerASCIIPrintf(viewer,"KSP solver for upper A00 in upper triangular factor \n");
230: PetscViewerASCIIPushTab(viewer);
231: if (jac->kspupper) {KSPView(jac->kspupper,viewer);}
232: else {PetscViewerASCIIPrintf(viewer," not yet available\n");}
233: PetscViewerASCIIPopTab(viewer);
234: }
235: PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");
236: PetscViewerASCIIPushTab(viewer);
237: if (jac->kspschur) {
238: KSPView(jac->kspschur,viewer);
239: } else {
240: PetscViewerASCIIPrintf(viewer," not yet available\n");
241: }
242: PetscViewerASCIIPopTab(viewer);
243: PetscViewerASCIIPopTab(viewer);
244: } else if (isdraw && jac->head) {
245: PetscDraw draw;
246: PetscReal x,y,w,wd,h;
247: PetscInt cnt = 2;
248: char str[32];
250: PetscViewerDrawGetDraw(viewer,0,&draw);
251: PetscDrawGetCurrentPoint(draw,&x,&y);
252: if (jac->kspupper != jac->head->ksp) cnt++;
253: w = 2*PetscMin(1.0 - x,x);
254: wd = w/(cnt + 1);
256: PetscSNPrintf(str,32,"Schur fact. %s",PCFieldSplitSchurFactTypes[jac->schurfactorization]);
257: PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);
258: y -= h;
259: if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_USER && !jac->schur_user) {
260: PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[PC_FIELDSPLIT_SCHUR_PRE_A11]);
261: } else {
262: PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[jac->schurpre]);
263: }
264: PetscDrawStringBoxed(draw,x+wd*(cnt-1)/2.0,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);
265: y -= h;
266: x = x - wd*(cnt-1)/2.0;
268: PetscDrawPushCurrentPoint(draw,x,y);
269: KSPView(jac->head->ksp,viewer);
270: PetscDrawPopCurrentPoint(draw);
271: if (jac->kspupper != jac->head->ksp) {
272: x += wd;
273: PetscDrawPushCurrentPoint(draw,x,y);
274: KSPView(jac->kspupper,viewer);
275: PetscDrawPopCurrentPoint(draw);
276: }
277: x += wd;
278: PetscDrawPushCurrentPoint(draw,x,y);
279: KSPView(jac->kspschur,viewer);
280: PetscDrawPopCurrentPoint(draw);
281: }
282: return(0);
283: }
285: static PetscErrorCode PCView_FieldSplit_GKB(PC pc,PetscViewer viewer)
286: {
287: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
288: PetscErrorCode ierr;
289: PetscBool iascii,isdraw;
290: PetscInt i,j;
291: PC_FieldSplitLink ilink = jac->head;
294: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
295: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
296: if (iascii) {
297: if (jac->bs > 0) {
298: PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);
299: } else {
300: PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D\n",PCCompositeTypes[jac->type],jac->nsplits);
301: }
302: if (pc->useAmat) {
303: PetscViewerASCIIPrintf(viewer," using Amat (not Pmat) as operator for blocks\n");
304: }
305: if (jac->diag_use_amat) {
306: PetscViewerASCIIPrintf(viewer," using Amat (not Pmat) as operator for diagonal blocks\n");
307: }
308: if (jac->offdiag_use_amat) {
309: PetscViewerASCIIPrintf(viewer," using Amat (not Pmat) as operator for off-diagonal blocks\n");
310: }
312: PetscViewerASCIIPrintf(viewer," Stopping tolerance=%.1e, delay in error estimate=%D, maximum iterations=%D\n",jac->gkbtol,jac->gkbdelay,jac->gkbmaxit);
313: PetscViewerASCIIPrintf(viewer," Solver info for H = A00 + nu*A01*A01' matrix:\n");
314: PetscViewerASCIIPushTab(viewer);
316: if (ilink->fields) {
317: PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",0);
318: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
319: for (j=0; j<ilink->nfields; j++) {
320: if (j > 0) {
321: PetscViewerASCIIPrintf(viewer,",");
322: }
323: PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);
324: }
325: PetscViewerASCIIPrintf(viewer,"\n");
326: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
327: } else {
328: PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",0);
329: }
330: KSPView(ilink->ksp,viewer);
332: PetscViewerASCIIPopTab(viewer);
333: }
335: if (isdraw) {
336: PetscDraw draw;
337: PetscReal x,y,w,wd;
339: PetscViewerDrawGetDraw(viewer,0,&draw);
340: PetscDrawGetCurrentPoint(draw,&x,&y);
341: w = 2*PetscMin(1.0 - x,x);
342: wd = w/(jac->nsplits + 1);
343: x = x - wd*(jac->nsplits-1)/2.0;
344: for (i=0; i<jac->nsplits; i++) {
345: PetscDrawPushCurrentPoint(draw,x,y);
346: KSPView(ilink->ksp,viewer);
347: PetscDrawPopCurrentPoint(draw);
348: x += wd;
349: ilink = ilink->next;
350: }
351: }
352: return(0);
353: }
356: /* Precondition: jac->bs is set to a meaningful value */
357: static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc)
358: {
360: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
361: PetscInt i,nfields,*ifields,nfields_col,*ifields_col;
362: PetscBool flg,flg_col;
363: char optionname[128],splitname[8],optionname_col[128];
366: PetscMalloc1(jac->bs,&ifields);
367: PetscMalloc1(jac->bs,&ifields_col);
368: for (i=0,flg=PETSC_TRUE;; i++) {
369: PetscSNPrintf(splitname,sizeof(splitname),"%D",i);
370: PetscSNPrintf(optionname,sizeof(optionname),"-pc_fieldsplit_%D_fields",i);
371: PetscSNPrintf(optionname_col,sizeof(optionname_col),"-pc_fieldsplit_%D_fields_col",i);
372: nfields = jac->bs;
373: nfields_col = jac->bs;
374: PetscOptionsGetIntArray(((PetscObject)pc)->options,((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);
375: PetscOptionsGetIntArray(((PetscObject)pc)->options,((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);
376: if (!flg) break;
377: else if (flg && !flg_col) {
378: if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
379: PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);
380: } else {
381: if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
382: if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match");
383: PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);
384: }
385: }
386: if (i > 0) {
387: /* Makes command-line setting of splits take precedence over setting them in code.
388: Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would
389: create new splits, which would probably not be what the user wanted. */
390: jac->splitdefined = PETSC_TRUE;
391: }
392: PetscFree(ifields);
393: PetscFree(ifields_col);
394: return(0);
395: }
397: static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
398: {
399: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
400: PetscErrorCode ierr;
401: PC_FieldSplitLink ilink = jac->head;
402: PetscBool fieldsplit_default = PETSC_FALSE,coupling = PETSC_FALSE;
403: PetscInt i;
406: /*
407: Kinda messy, but at least this now uses DMCreateFieldDecomposition().
408: Should probably be rewritten.
409: */
410: if (!ilink) {
411: PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_coupling",&coupling,NULL);
412: if (pc->dm && jac->dm_splits && !jac->detect && !coupling) {
413: PetscInt numFields, f, i, j;
414: char **fieldNames;
415: IS *fields;
416: DM *dms;
417: DM subdm[128];
418: PetscBool flg;
420: DMCreateFieldDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms);
421: /* Allow the user to prescribe the splits */
422: for (i = 0, flg = PETSC_TRUE;; i++) {
423: PetscInt ifields[128];
424: IS compField;
425: char optionname[128], splitname[8];
426: PetscInt nfields = numFields;
428: PetscSNPrintf(optionname, sizeof(optionname), "-pc_fieldsplit_%D_fields", i);
429: PetscOptionsGetIntArray(((PetscObject)pc)->options,((PetscObject)pc)->prefix, optionname, ifields, &nfields, &flg);
430: if (!flg) break;
431: if (numFields > 128) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot currently support %d > 128 fields", numFields);
432: DMCreateSubDM(pc->dm, nfields, ifields, &compField, &subdm[i]);
433: if (nfields == 1) {
434: PCFieldSplitSetIS(pc, fieldNames[ifields[0]], compField);
435: } else {
436: PetscSNPrintf(splitname, sizeof(splitname), "%D", i);
437: PCFieldSplitSetIS(pc, splitname, compField);
438: }
439: ISDestroy(&compField);
440: for (j = 0; j < nfields; ++j) {
441: f = ifields[j];
442: PetscFree(fieldNames[f]);
443: ISDestroy(&fields[f]);
444: }
445: }
446: if (i == 0) {
447: for (f = 0; f < numFields; ++f) {
448: PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);
449: PetscFree(fieldNames[f]);
450: ISDestroy(&fields[f]);
451: }
452: } else {
453: for (j=0; j<numFields; j++) {
454: DMDestroy(dms+j);
455: }
456: PetscFree(dms);
457: PetscMalloc1(i, &dms);
458: for (j = 0; j < i; ++j) dms[j] = subdm[j];
459: }
460: PetscFree(fieldNames);
461: PetscFree(fields);
462: if (dms) {
463: PetscInfo(pc, "Setting up physics based fieldsplit preconditioner using the embedded DM\n");
464: for (ilink = jac->head, i = 0; ilink; ilink = ilink->next, ++i) {
465: const char *prefix;
466: PetscObjectGetOptionsPrefix((PetscObject)(ilink->ksp),&prefix);
467: PetscObjectSetOptionsPrefix((PetscObject)(dms[i]), prefix);
468: KSPSetDM(ilink->ksp, dms[i]);
469: KSPSetDMActive(ilink->ksp, PETSC_FALSE);
470: {
471: PetscErrorCode (*func)(KSP,Mat,Mat,void*);
472: void *ctx;
474: DMKSPGetComputeOperators(pc->dm, &func, &ctx);
475: DMKSPSetComputeOperators(dms[i], func, ctx);
476: }
477: PetscObjectIncrementTabLevel((PetscObject)dms[i],(PetscObject)ilink->ksp,0);
478: DMDestroy(&dms[i]);
479: }
480: PetscFree(dms);
481: }
482: } else {
483: if (jac->bs <= 0) {
484: if (pc->pmat) {
485: MatGetBlockSize(pc->pmat,&jac->bs);
486: } else jac->bs = 1;
487: }
489: if (jac->detect) {
490: IS zerodiags,rest;
491: PetscInt nmin,nmax;
493: MatGetOwnershipRange(pc->mat,&nmin,&nmax);
494: MatFindZeroDiagonals(pc->mat,&zerodiags);
495: ISComplement(zerodiags,nmin,nmax,&rest);
496: PCFieldSplitSetIS(pc,"0",rest);
497: PCFieldSplitSetIS(pc,"1",zerodiags);
498: ISDestroy(&zerodiags);
499: ISDestroy(&rest);
500: } else if (coupling) {
501: IS coupling,rest;
502: PetscInt nmin,nmax;
504: MatGetOwnershipRange(pc->mat,&nmin,&nmax);
505: MatFindOffBlockDiagonalEntries(pc->mat,&coupling);
506: ISCreateStride(PetscObjectComm((PetscObject)pc->mat),nmax-nmin,nmin,1,&rest);
507: ISSetIdentity(rest);
508: PCFieldSplitSetIS(pc,"0",rest);
509: PCFieldSplitSetIS(pc,"1",coupling);
510: ISDestroy(&coupling);
511: ISDestroy(&rest);
512: } else {
513: PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&fieldsplit_default,NULL);
514: if (!fieldsplit_default) {
515: /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit()
516: then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
517: PCFieldSplitSetRuntimeSplits_Private(pc);
518: if (jac->splitdefined) {PetscInfo(pc,"Splits defined using the options database\n");}
519: }
520: if ((fieldsplit_default || !jac->splitdefined) && !jac->isrestrict) {
521: Mat M = pc->pmat;
522: PetscBool isnest;
524: PetscInfo(pc,"Using default splitting of fields\n");
525: PetscObjectTypeCompare((PetscObject)pc->pmat,MATNEST,&isnest);
526: if (!isnest) {
527: M = pc->mat;
528: PetscObjectTypeCompare((PetscObject)pc->mat,MATNEST,&isnest);
529: }
530: if (isnest) {
531: IS *fields;
532: PetscInt nf;
534: MatNestGetSize(M,&nf,NULL);
535: PetscMalloc1(nf,&fields);
536: MatNestGetISs(M,fields,NULL);
537: for (i=0;i<nf;i++) {
538: PCFieldSplitSetIS(pc,NULL,fields[i]);
539: }
540: PetscFree(fields);
541: } else {
542: for (i=0; i<jac->bs; i++) {
543: char splitname[8];
544: PetscSNPrintf(splitname,sizeof(splitname),"%D",i);
545: PCFieldSplitSetFields(pc,splitname,1,&i,&i);
546: }
547: jac->defaultsplit = PETSC_TRUE;
548: }
549: }
550: }
551: }
552: } else if (jac->nsplits == 1) {
553: if (ilink->is) {
554: IS is2;
555: PetscInt nmin,nmax;
557: MatGetOwnershipRange(pc->mat,&nmin,&nmax);
558: ISComplement(ilink->is,nmin,nmax,&is2);
559: PCFieldSplitSetIS(pc,"1",is2);
560: ISDestroy(&is2);
561: } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
562: }
564: if (jac->nsplits < 2) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits);
565: return(0);
566: }
568: static PetscErrorCode MatGolubKahanComputeExplicitOperator(Mat A,Mat B,Mat C,Mat *H,PetscReal gkbnu)
569: {
570: PetscErrorCode ierr;
571: Mat BT,T;
572: PetscReal nrmT,nrmB;
575: MatHermitianTranspose(C,MAT_INITIAL_MATRIX,&T); /* Test if augmented matrix is symmetric */
576: MatAXPY(T,-1.0,B,DIFFERENT_NONZERO_PATTERN);
577: MatNorm(T,NORM_1,&nrmT);
578: MatNorm(B,NORM_1,&nrmB);
579: if (nrmB > 0) {
580: if (nrmT/nrmB >= PETSC_SMALL) {
581: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix is not symmetric/hermitian, GKB is not applicable.");
582: }
583: }
584: /* Compute augmented Lagrangian matrix H = A00 + nu*A01*A01'. This corresponds to */
585: /* setting N := 1/nu*I in [Ar13]. */
586: MatHermitianTranspose(B,MAT_INITIAL_MATRIX,&BT);
587: MatMatMult(B,BT,MAT_INITIAL_MATRIX,PETSC_DEFAULT,H); /* H = A01*A01' */
588: MatAYPX(*H,gkbnu,A,DIFFERENT_NONZERO_PATTERN); /* H = A00 + nu*A01*A01' */
590: MatDestroy(&BT);
591: MatDestroy(&T);
592: return(0);
593: }
595: PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(PetscOptions,const char pre[], const char name[],const char *value[],PetscBool *flg);
597: static PetscErrorCode PCSetUp_FieldSplit(PC pc)
598: {
599: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
600: PetscErrorCode ierr;
601: PC_FieldSplitLink ilink;
602: PetscInt i,nsplit;
603: PetscBool sorted, sorted_col;
606: pc->failedreason = PC_NOERROR;
607: PCFieldSplitSetDefaults(pc);
608: nsplit = jac->nsplits;
609: ilink = jac->head;
611: /* get the matrices for each split */
612: if (!jac->issetup) {
613: PetscInt rstart,rend,nslots,bs;
615: jac->issetup = PETSC_TRUE;
617: /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
618: if (jac->defaultsplit || !ilink->is) {
619: if (jac->bs <= 0) jac->bs = nsplit;
620: }
621: bs = jac->bs;
622: MatGetOwnershipRange(pc->pmat,&rstart,&rend);
623: nslots = (rend - rstart)/bs;
624: for (i=0; i<nsplit; i++) {
625: if (jac->defaultsplit) {
626: ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+i,nsplit,&ilink->is);
627: ISDuplicate(ilink->is,&ilink->is_col);
628: } else if (!ilink->is) {
629: if (ilink->nfields > 1) {
630: PetscInt *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col;
631: PetscMalloc1(ilink->nfields*nslots,&ii);
632: PetscMalloc1(ilink->nfields*nslots,&jj);
633: for (j=0; j<nslots; j++) {
634: for (k=0; k<nfields; k++) {
635: ii[nfields*j + k] = rstart + bs*j + fields[k];
636: jj[nfields*j + k] = rstart + bs*j + fields_col[k];
637: }
638: }
639: ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);
640: ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);
641: ISSetBlockSize(ilink->is, nfields);
642: ISSetBlockSize(ilink->is_col, nfields);
643: } else {
644: ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields[0],bs,&ilink->is);
645: ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);
646: }
647: }
648: ISSorted(ilink->is,&sorted);
649: if (ilink->is_col) { ISSorted(ilink->is_col,&sorted_col); }
650: if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
651: ilink = ilink->next;
652: }
653: }
655: ilink = jac->head;
656: if (!jac->pmat) {
657: Vec xtmp;
659: MatCreateVecs(pc->pmat,&xtmp,NULL);
660: PetscMalloc1(nsplit,&jac->pmat);
661: PetscMalloc2(nsplit,&jac->x,nsplit,&jac->y);
662: for (i=0; i<nsplit; i++) {
663: MatNullSpace sp;
665: /* Check for preconditioning matrix attached to IS */
666: PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &jac->pmat[i]);
667: if (jac->pmat[i]) {
668: PetscObjectReference((PetscObject) jac->pmat[i]);
669: if (jac->type == PC_COMPOSITE_SCHUR) {
670: jac->schur_user = jac->pmat[i];
672: PetscObjectReference((PetscObject) jac->schur_user);
673: }
674: } else {
675: const char *prefix;
676: MatCreateSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);
677: KSPGetOptionsPrefix(ilink->ksp,&prefix);
678: MatSetOptionsPrefix(jac->pmat[i],prefix);
679: MatViewFromOptions(jac->pmat[i],NULL,"-mat_view");
680: }
681: /* create work vectors for each split */
682: MatCreateVecs(jac->pmat[i],&jac->x[i],&jac->y[i]);
683: ilink->x = jac->x[i]; ilink->y = jac->y[i]; ilink->z = NULL;
684: /* compute scatter contexts needed by multiplicative versions and non-default splits */
685: VecScatterCreate(xtmp,ilink->is,jac->x[i],NULL,&ilink->sctx);
686: PetscObjectQuery((PetscObject) ilink->is, "nearnullspace", (PetscObject*) &sp);
687: if (sp) {
688: MatSetNearNullSpace(jac->pmat[i], sp);
689: }
690: ilink = ilink->next;
691: }
692: VecDestroy(&xtmp);
693: } else {
694: MatReuse scall;
695: if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
696: for (i=0; i<nsplit; i++) {
697: MatDestroy(&jac->pmat[i]);
698: }
699: scall = MAT_INITIAL_MATRIX;
700: } else scall = MAT_REUSE_MATRIX;
702: for (i=0; i<nsplit; i++) {
703: Mat pmat;
705: /* Check for preconditioning matrix attached to IS */
706: PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &pmat);
707: if (!pmat) {
708: MatCreateSubMatrix(pc->pmat,ilink->is,ilink->is_col,scall,&jac->pmat[i]);
709: }
710: ilink = ilink->next;
711: }
712: }
713: if (jac->diag_use_amat) {
714: ilink = jac->head;
715: if (!jac->mat) {
716: PetscMalloc1(nsplit,&jac->mat);
717: for (i=0; i<nsplit; i++) {
718: MatCreateSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);
719: ilink = ilink->next;
720: }
721: } else {
722: MatReuse scall;
723: if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
724: for (i=0; i<nsplit; i++) {
725: MatDestroy(&jac->mat[i]);
726: }
727: scall = MAT_INITIAL_MATRIX;
728: } else scall = MAT_REUSE_MATRIX;
730: for (i=0; i<nsplit; i++) {
731: MatCreateSubMatrix(pc->mat,ilink->is,ilink->is_col,scall,&jac->mat[i]);
732: ilink = ilink->next;
733: }
734: }
735: } else {
736: jac->mat = jac->pmat;
737: }
739: /* Check for null space attached to IS */
740: ilink = jac->head;
741: for (i=0; i<nsplit; i++) {
742: MatNullSpace sp;
744: PetscObjectQuery((PetscObject) ilink->is, "nullspace", (PetscObject*) &sp);
745: if (sp) {
746: MatSetNullSpace(jac->mat[i], sp);
747: }
748: ilink = ilink->next;
749: }
751: if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR && jac->type != PC_COMPOSITE_GKB) {
752: /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
753: /* FIXME: Can/should we reuse jac->mat whenever (jac->diag_use_amat) is true? */
754: ilink = jac->head;
755: if (nsplit == 2 && jac->type == PC_COMPOSITE_MULTIPLICATIVE) {
756: /* special case need where Afield[0] is not needed and only certain columns of Afield[1] are needed since update is only on those rows of the solution */
757: if (!jac->Afield) {
758: PetscCalloc1(nsplit,&jac->Afield);
759: if (jac->offdiag_use_amat) {
760: MatCreateSubMatrix(pc->mat,ilink->next->is,ilink->is,MAT_INITIAL_MATRIX,&jac->Afield[1]);
761: } else {
762: MatCreateSubMatrix(pc->pmat,ilink->next->is,ilink->is,MAT_INITIAL_MATRIX,&jac->Afield[1]);
763: }
764: } else {
765: MatReuse scall;
767: if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
768: MatDestroy(&jac->Afield[1]);
769: scall = MAT_INITIAL_MATRIX;
770: } else scall = MAT_REUSE_MATRIX;
772: if (jac->offdiag_use_amat) {
773: MatCreateSubMatrix(pc->mat,ilink->next->is,ilink->is,scall,&jac->Afield[1]);
774: } else {
775: MatCreateSubMatrix(pc->pmat,ilink->next->is,ilink->is,scall,&jac->Afield[1]);
776: }
777: }
778: } else {
779: if (!jac->Afield) {
780: PetscMalloc1(nsplit,&jac->Afield);
781: for (i=0; i<nsplit; i++) {
782: if (jac->offdiag_use_amat) {
783: MatCreateSubMatrix(pc->mat,ilink->is,NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);
784: } else {
785: MatCreateSubMatrix(pc->pmat,ilink->is,NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);
786: }
787: ilink = ilink->next;
788: }
789: } else {
790: MatReuse scall;
791: if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
792: for (i=0; i<nsplit; i++) {
793: MatDestroy(&jac->Afield[i]);
794: }
795: scall = MAT_INITIAL_MATRIX;
796: } else scall = MAT_REUSE_MATRIX;
798: for (i=0; i<nsplit; i++) {
799: if (jac->offdiag_use_amat) {
800: MatCreateSubMatrix(pc->mat,ilink->is,NULL,scall,&jac->Afield[i]);
801: } else {
802: MatCreateSubMatrix(pc->pmat,ilink->is,NULL,scall,&jac->Afield[i]);
803: }
804: ilink = ilink->next;
805: }
806: }
807: }
808: }
810: if (jac->type == PC_COMPOSITE_SCHUR) {
811: IS ccis;
812: PetscBool isspd;
813: PetscInt rstart,rend;
814: char lscname[256];
815: PetscObject LSC_L;
817: if (nsplit != 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
819: /* If pc->mat is SPD, don't scale by -1 the Schur complement */
820: if (jac->schurscale == (PetscScalar)-1.0) {
821: MatGetOption(pc->pmat,MAT_SPD,&isspd);
822: jac->schurscale = (isspd == PETSC_TRUE) ? 1.0 : -1.0;
823: }
825: /* When extracting off-diagonal submatrices, we take complements from this range */
826: MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);
828: if (jac->schur) {
829: KSP kspA = jac->head->ksp, kspInner = NULL, kspUpper = jac->kspupper;
830: MatReuse scall;
832: if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
833: scall = MAT_INITIAL_MATRIX;
834: MatDestroy(&jac->B);
835: MatDestroy(&jac->C);
836: } else scall = MAT_REUSE_MATRIX;
838: MatSchurComplementGetKSP(jac->schur, &kspInner);
839: ilink = jac->head;
840: ISComplement(ilink->is_col,rstart,rend,&ccis);
841: if (jac->offdiag_use_amat) {
842: MatCreateSubMatrix(pc->mat,ilink->is,ccis,scall,&jac->B);
843: } else {
844: MatCreateSubMatrix(pc->pmat,ilink->is,ccis,scall,&jac->B);
845: }
846: ISDestroy(&ccis);
847: ilink = ilink->next;
848: ISComplement(ilink->is_col,rstart,rend,&ccis);
849: if (jac->offdiag_use_amat) {
850: MatCreateSubMatrix(pc->mat,ilink->is,ccis,scall,&jac->C);
851: } else {
852: MatCreateSubMatrix(pc->pmat,ilink->is,ccis,scall,&jac->C);
853: }
854: ISDestroy(&ccis);
855: MatSchurComplementUpdateSubMatrices(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);
856: if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELFP) {
857: MatDestroy(&jac->schurp);
858: MatSchurComplementGetPmat(jac->schur,MAT_INITIAL_MATRIX,&jac->schurp);
859: }
860: if (kspA != kspInner) {
861: KSPSetOperators(kspA,jac->mat[0],jac->pmat[0]);
862: }
863: if (kspUpper != kspA) {
864: KSPSetOperators(kspUpper,jac->mat[0],jac->pmat[0]);
865: }
866: KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac));
867: } else {
868: const char *Dprefix;
869: char schurprefix[256], schurmatprefix[256];
870: char schurtestoption[256];
871: MatNullSpace sp;
872: PetscBool flg;
873: KSP kspt;
875: /* extract the A01 and A10 matrices */
876: ilink = jac->head;
877: ISComplement(ilink->is_col,rstart,rend,&ccis);
878: if (jac->offdiag_use_amat) {
879: MatCreateSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);
880: } else {
881: MatCreateSubMatrix(pc->pmat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);
882: }
883: ISDestroy(&ccis);
884: ilink = ilink->next;
885: ISComplement(ilink->is_col,rstart,rend,&ccis);
886: if (jac->offdiag_use_amat) {
887: MatCreateSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);
888: } else {
889: MatCreateSubMatrix(pc->pmat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);
890: }
891: ISDestroy(&ccis);
893: /* Use mat[0] (diagonal block of Amat) preconditioned by pmat[0] to define Schur complement */
894: MatCreate(((PetscObject)jac->mat[0])->comm,&jac->schur);
895: MatSetType(jac->schur,MATSCHURCOMPLEMENT);
896: MatSchurComplementSetSubMatrices(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);
897: PetscSNPrintf(schurmatprefix, sizeof(schurmatprefix), "%sfieldsplit_%s_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);
898: MatSetOptionsPrefix(jac->schur,schurmatprefix);
899: MatSchurComplementGetKSP(jac->schur,&kspt);
900: KSPSetOptionsPrefix(kspt,schurmatprefix);
902: /* Note: this is not true in general */
903: MatGetNullSpace(jac->mat[1], &sp);
904: if (sp) {
905: MatSetNullSpace(jac->schur, sp);
906: }
908: PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_inner_", ilink->splitname);
909: PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->options,((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);
910: if (flg) {
911: DM dmInner;
912: KSP kspInner;
913: PC pcInner;
915: MatSchurComplementGetKSP(jac->schur, &kspInner);
916: KSPReset(kspInner);
917: KSPSetOperators(kspInner,jac->mat[0],jac->pmat[0]);
918: PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_inner_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);
919: /* Indent this deeper to emphasize the "inner" nature of this solver. */
920: PetscObjectIncrementTabLevel((PetscObject)kspInner, (PetscObject) pc, 2);
921: PetscObjectIncrementTabLevel((PetscObject)kspInner->pc, (PetscObject) pc, 2);
922: KSPSetOptionsPrefix(kspInner, schurprefix);
924: /* Set DM for new solver */
925: KSPGetDM(jac->head->ksp, &dmInner);
926: KSPSetDM(kspInner, dmInner);
927: KSPSetDMActive(kspInner, PETSC_FALSE);
929: /* Defaults to PCKSP as preconditioner */
930: KSPGetPC(kspInner, &pcInner);
931: PCSetType(pcInner, PCKSP);
932: PCKSPSetKSP(pcInner, jac->head->ksp);
933: } else {
934: /* Use the outer solver for the inner solve, but revert the KSPPREONLY from PCFieldSplitSetFields_FieldSplit or
935: * PCFieldSplitSetIS_FieldSplit. We don't want KSPPREONLY because it makes the Schur complement inexact,
936: * preventing Schur complement reduction to be an accurate solve. Usually when an iterative solver is used for
937: * S = D - C A_inner^{-1} B, we expect S to be defined using an accurate definition of A_inner^{-1}, so we make
938: * GMRES the default. Note that it is also common to use PREONLY for S, in which case S may not be used
939: * directly, and the user is responsible for setting an inexact method for fieldsplit's A^{-1}. */
940: KSPSetType(jac->head->ksp,KSPGMRES);
941: MatSchurComplementSetKSP(jac->schur,jac->head->ksp);
942: }
943: KSPSetOperators(jac->head->ksp,jac->mat[0],jac->pmat[0]);
944: KSPSetFromOptions(jac->head->ksp);
945: MatSetFromOptions(jac->schur);
947: PetscObjectTypeCompare((PetscObject)jac->schur, MATSCHURCOMPLEMENT, &flg);
948: if (flg) { /* Need to do this otherwise PCSetUp_KSP will overwrite the amat of jac->head->ksp */
949: KSP kspInner;
950: PC pcInner;
952: MatSchurComplementGetKSP(jac->schur, &kspInner);
953: KSPGetPC(kspInner, &pcInner);
954: PetscObjectTypeCompare((PetscObject)pcInner, PCKSP, &flg);
955: if (flg) {
956: KSP ksp;
958: PCKSPGetKSP(pcInner, &ksp);
959: if (ksp == jac->head->ksp) {
960: PCSetUseAmat(pcInner, PETSC_TRUE);
961: }
962: }
963: }
964: PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_upper_", ilink->splitname);
965: PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->options,((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);
966: if (flg) {
967: DM dmInner;
969: PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_upper_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);
970: KSPCreate(PetscObjectComm((PetscObject)pc), &jac->kspupper);
971: KSPSetErrorIfNotConverged(jac->kspupper,pc->erroriffailure);
972: KSPSetOptionsPrefix(jac->kspupper, schurprefix);
973: PetscObjectIncrementTabLevel((PetscObject)jac->kspupper, (PetscObject) pc, 1);
974: PetscObjectIncrementTabLevel((PetscObject)jac->kspupper->pc, (PetscObject) pc, 1);
975: KSPGetDM(jac->head->ksp, &dmInner);
976: KSPSetDM(jac->kspupper, dmInner);
977: KSPSetDMActive(jac->kspupper, PETSC_FALSE);
978: KSPSetFromOptions(jac->kspupper);
979: KSPSetOperators(jac->kspupper,jac->mat[0],jac->pmat[0]);
980: VecDuplicate(jac->head->x, &jac->head->z);
981: } else {
982: jac->kspupper = jac->head->ksp;
983: PetscObjectReference((PetscObject) jac->head->ksp);
984: }
986: if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELFP) {
987: MatSchurComplementGetPmat(jac->schur,MAT_INITIAL_MATRIX,&jac->schurp);
988: }
989: KSPCreate(PetscObjectComm((PetscObject)pc),&jac->kspschur);
990: KSPSetErrorIfNotConverged(jac->kspschur,pc->erroriffailure);
991: PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);
992: PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);
993: if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
994: PC pcschur;
995: KSPGetPC(jac->kspschur,&pcschur);
996: PCSetType(pcschur,PCNONE);
997: /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
998: } else if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_FULL) {
999: MatSchurComplementComputeExplicitOperator(jac->schur, &jac->schur_user);
1000: }
1001: KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac));
1002: KSPGetOptionsPrefix(jac->head->next->ksp, &Dprefix);
1003: KSPSetOptionsPrefix(jac->kspschur, Dprefix);
1004: /* propagate DM */
1005: {
1006: DM sdm;
1007: KSPGetDM(jac->head->next->ksp, &sdm);
1008: if (sdm) {
1009: KSPSetDM(jac->kspschur, sdm);
1010: KSPSetDMActive(jac->kspschur, PETSC_FALSE);
1011: }
1012: }
1013: /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
1014: /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */
1015: KSPSetFromOptions(jac->kspschur);
1016: }
1017: MatAssemblyBegin(jac->schur,MAT_FINAL_ASSEMBLY);
1018: MatAssemblyEnd(jac->schur,MAT_FINAL_ASSEMBLY);
1020: /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */
1021: PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_L",ilink->splitname);
1022: PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);
1023: if (!LSC_L) {PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);}
1024: if (LSC_L) {PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);}
1025: PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_Lp",ilink->splitname);
1026: PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);
1027: if (!LSC_L) {PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);}
1028: if (LSC_L) {PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);}
1029: } else if (jac->type == PC_COMPOSITE_GKB) {
1030: IS ccis;
1031: PetscInt rstart,rend;
1033: if (nsplit != 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"To use GKB preconditioner you must have exactly 2 fields");
1035: ilink = jac->head;
1037: /* When extracting off-diagonal submatrices, we take complements from this range */
1038: MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);
1040: ISComplement(ilink->is_col,rstart,rend,&ccis);
1041: if (jac->offdiag_use_amat) {
1042: MatCreateSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);
1043: } else {
1044: MatCreateSubMatrix(pc->pmat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);
1045: }
1046: ISDestroy(&ccis);
1047: /* Create work vectors for GKB algorithm */
1048: VecDuplicate(ilink->x,&jac->u);
1049: VecDuplicate(ilink->x,&jac->Hu);
1050: VecDuplicate(ilink->x,&jac->w2);
1051: ilink = ilink->next;
1052: ISComplement(ilink->is_col,rstart,rend,&ccis);
1053: if (jac->offdiag_use_amat) {
1054: MatCreateSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);
1055: } else {
1056: MatCreateSubMatrix(pc->pmat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);
1057: }
1058: ISDestroy(&ccis);
1059: /* Create work vectors for GKB algorithm */
1060: VecDuplicate(ilink->x,&jac->v);
1061: VecDuplicate(ilink->x,&jac->d);
1062: VecDuplicate(ilink->x,&jac->w1);
1063: MatGolubKahanComputeExplicitOperator(jac->mat[0],jac->B,jac->C,&jac->H,jac->gkbnu);
1064: PetscCalloc1(jac->gkbdelay,&jac->vecz);
1066: ilink = jac->head;
1067: KSPSetOperators(ilink->ksp,jac->H,jac->H);
1068: if (!jac->suboptionsset) {KSPSetFromOptions(ilink->ksp);}
1069: /* Create gkb_monitor context */
1070: if (jac->gkbmonitor) {
1071: PetscInt tablevel;
1072: PetscViewerCreate(PETSC_COMM_WORLD,&jac->gkbviewer);
1073: PetscViewerSetType(jac->gkbviewer,PETSCVIEWERASCII);
1074: PetscObjectGetTabLevel((PetscObject)ilink->ksp,&tablevel);
1075: PetscViewerASCIISetTab(jac->gkbviewer,tablevel);
1076: PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)ilink->ksp,1);
1077: }
1078: } else {
1079: /* set up the individual splits' PCs */
1080: i = 0;
1081: ilink = jac->head;
1082: while (ilink) {
1083: KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i]);
1084: /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
1085: if (!jac->suboptionsset) {KSPSetFromOptions(ilink->ksp);}
1086: i++;
1087: ilink = ilink->next;
1088: }
1089: }
1091: jac->suboptionsset = PETSC_TRUE;
1092: return(0);
1093: }
1095: #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
1096: (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
1097: VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
1098: PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL) ||\
1099: KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
1100: KSPCheckSolve(ilink->ksp,pc,ilink->y) || \
1101: PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL) ||\
1102: VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
1103: VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
1105: static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
1106: {
1107: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1108: PetscErrorCode ierr;
1109: PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
1110: KSP kspA = ilinkA->ksp, kspLower = kspA, kspUpper = jac->kspupper;
1113: switch (jac->schurfactorization) {
1114: case PC_FIELDSPLIT_SCHUR_FACT_DIAG:
1115: /* [A00 0; 0 -S], positive definite, suitable for MINRES */
1116: VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
1117: VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);
1118: VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
1119: PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
1120: KSPSolve(kspA,ilinkA->x,ilinkA->y);
1121: KSPCheckSolve(kspA,pc,ilinkA->y);
1122: PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
1123: VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
1124: VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);
1125: PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);
1126: KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);
1127: KSPCheckSolve(jac->kspschur,pc,ilinkD->y);
1128: PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);
1129: VecScale(ilinkD->y,jac->schurscale);
1130: VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
1131: VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
1132: VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
1133: break;
1134: case PC_FIELDSPLIT_SCHUR_FACT_LOWER:
1135: /* [A00 0; A10 S], suitable for left preconditioning */
1136: VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
1137: VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
1138: PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
1139: KSPSolve(kspA,ilinkA->x,ilinkA->y);
1140: KSPCheckSolve(kspA,pc,ilinkA->y);
1141: PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
1142: MatMult(jac->C,ilinkA->y,ilinkD->x);
1143: VecScale(ilinkD->x,-1.);
1144: VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);
1145: VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
1146: VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);
1147: PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);
1148: KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);
1149: KSPCheckSolve(jac->kspschur,pc,ilinkD->y);
1150: PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);
1151: VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
1152: VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
1153: VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
1154: break;
1155: case PC_FIELDSPLIT_SCHUR_FACT_UPPER:
1156: /* [A00 A01; 0 S], suitable for right preconditioning */
1157: VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);
1158: VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);
1159: PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);
1160: KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);
1161: KSPCheckSolve(jac->kspschur,pc,ilinkD->y);
1162: PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL); MatMult(jac->B,ilinkD->y,ilinkA->x);
1163: VecScale(ilinkA->x,-1.);
1164: VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);
1165: VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
1166: VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);
1167: PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
1168: KSPSolve(kspA,ilinkA->x,ilinkA->y);
1169: KSPCheckSolve(kspA,pc,ilinkA->y);
1170: PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
1171: VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
1172: VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
1173: VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
1174: break;
1175: case PC_FIELDSPLIT_SCHUR_FACT_FULL:
1176: /* [1 0; A10 A00^{-1} 1] [A00 0; 0 S] [1 A00^{-1}A01; 0 1] */
1177: VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
1178: VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
1179: PetscLogEventBegin(KSP_Solve_FS_L,kspLower,ilinkA->x,ilinkA->y,NULL);
1180: KSPSolve(kspLower,ilinkA->x,ilinkA->y);
1181: KSPCheckSolve(kspLower,pc,ilinkA->y);
1182: PetscLogEventEnd(KSP_Solve_FS_L,kspLower,ilinkA->x,ilinkA->y,NULL);
1183: MatMult(jac->C,ilinkA->y,ilinkD->x);
1184: VecScale(ilinkD->x,-1.0);
1185: VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);
1186: VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);
1188: PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);
1189: KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);
1190: KSPCheckSolve(jac->kspschur,pc,ilinkD->y);
1191: PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);
1192: VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
1194: if (kspUpper == kspA) {
1195: MatMult(jac->B,ilinkD->y,ilinkA->y);
1196: VecAXPY(ilinkA->x,-1.0,ilinkA->y);
1197: PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
1198: KSPSolve(kspA,ilinkA->x,ilinkA->y);
1199: KSPCheckSolve(kspA,pc,ilinkA->y);
1200: PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
1201: } else {
1202: PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);
1203: KSPSolve(kspA,ilinkA->x,ilinkA->y);
1204: KSPCheckSolve(kspA,pc,ilinkA->y);
1205: MatMult(jac->B,ilinkD->y,ilinkA->x);
1206: PetscLogEventBegin(KSP_Solve_FS_U,kspUpper,ilinkA->x,ilinkA->z,NULL);
1207: KSPSolve(kspUpper,ilinkA->x,ilinkA->z);
1208: KSPCheckSolve(kspUpper,pc,ilinkA->z);
1209: PetscLogEventEnd(KSP_Solve_FS_U,kspUpper,ilinkA->x,ilinkA->z,NULL);
1210: VecAXPY(ilinkA->y,-1.0,ilinkA->z);
1211: }
1212: VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
1213: VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
1214: VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
1215: }
1216: return(0);
1217: }
1219: static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
1220: {
1221: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1222: PetscErrorCode ierr;
1223: PC_FieldSplitLink ilink = jac->head;
1224: PetscInt cnt,bs;
1227: if (jac->type == PC_COMPOSITE_ADDITIVE) {
1228: if (jac->defaultsplit) {
1229: VecGetBlockSize(x,&bs);
1230: if (jac->bs > 0 && bs != jac->bs) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Blocksize of x vector %D does not match fieldsplit blocksize %D",bs,jac->bs);
1231: VecGetBlockSize(y,&bs);
1232: if (jac->bs > 0 && bs != jac->bs) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Blocksize of y vector %D does not match fieldsplit blocksize %D",bs,jac->bs);
1233: VecStrideGatherAll(x,jac->x,INSERT_VALUES);
1234: while (ilink) {
1235: PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1236: KSPSolve(ilink->ksp,ilink->x,ilink->y);
1237: KSPCheckSolve(ilink->ksp,pc,ilink->y);
1238: PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1239: ilink = ilink->next;
1240: }
1241: VecStrideScatterAll(jac->y,y,INSERT_VALUES);
1242: } else {
1243: VecSet(y,0.0);
1244: while (ilink) {
1245: FieldSplitSplitSolveAdd(ilink,x,y);
1246: ilink = ilink->next;
1247: }
1248: }
1249: } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE && jac->nsplits == 2) {
1250: VecSet(y,0.0);
1251: /* solve on first block for first block variables */
1252: VecScatterBegin(ilink->sctx,x,ilink->x,INSERT_VALUES,SCATTER_FORWARD);
1253: VecScatterEnd(ilink->sctx,x,ilink->x,INSERT_VALUES,SCATTER_FORWARD);
1254: PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1255: KSPSolve(ilink->ksp,ilink->x,ilink->y);
1256: KSPCheckSolve(ilink->ksp,pc,ilink->y);
1257: PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1258: VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
1259: VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
1261: /* compute the residual only onto second block variables using first block variables */
1262: MatMult(jac->Afield[1],ilink->y,ilink->next->x);
1263: ilink = ilink->next;
1264: VecScale(ilink->x,-1.0);
1265: VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);
1266: VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);
1268: /* solve on second block variables */
1269: PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1270: KSPSolve(ilink->ksp,ilink->x,ilink->y);
1271: KSPCheckSolve(ilink->ksp,pc,ilink->y);
1272: PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1273: VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
1274: VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
1275: } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
1276: if (!jac->w1) {
1277: VecDuplicate(x,&jac->w1);
1278: VecDuplicate(x,&jac->w2);
1279: }
1280: VecSet(y,0.0);
1281: FieldSplitSplitSolveAdd(ilink,x,y);
1282: cnt = 1;
1283: while (ilink->next) {
1284: ilink = ilink->next;
1285: /* compute the residual only over the part of the vector needed */
1286: MatMult(jac->Afield[cnt++],y,ilink->x);
1287: VecScale(ilink->x,-1.0);
1288: VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);
1289: VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);
1290: PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1291: KSPSolve(ilink->ksp,ilink->x,ilink->y);
1292: KSPCheckSolve(ilink->ksp,pc,ilink->y);
1293: PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1294: VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
1295: VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
1296: }
1297: if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
1298: cnt -= 2;
1299: while (ilink->previous) {
1300: ilink = ilink->previous;
1301: /* compute the residual only over the part of the vector needed */
1302: MatMult(jac->Afield[cnt--],y,ilink->x);
1303: VecScale(ilink->x,-1.0);
1304: VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);
1305: VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);
1306: PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1307: KSPSolve(ilink->ksp,ilink->x,ilink->y);
1308: KSPCheckSolve(ilink->ksp,pc,ilink->y);
1309: PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1310: VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
1311: VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);
1312: }
1313: }
1314: } else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
1315: return(0);
1316: }
1319: static PetscErrorCode PCApply_FieldSplit_GKB(PC pc,Vec x,Vec y)
1320: {
1321: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1322: PetscErrorCode ierr;
1323: PC_FieldSplitLink ilinkA = jac->head,ilinkD = ilinkA->next;
1324: KSP ksp = ilinkA->ksp;
1325: Vec u,v,Hu,d,work1,work2;
1326: PetscScalar alpha,z,nrmz2,*vecz;
1327: PetscReal lowbnd,nu,beta;
1328: PetscInt j,iterGKB;
1331: VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
1332: VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);
1333: VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);
1334: VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);
1336: u = jac->u;
1337: v = jac->v;
1338: Hu = jac->Hu;
1339: d = jac->d;
1340: work1 = jac->w1;
1341: work2 = jac->w2;
1342: vecz = jac->vecz;
1344: /* Change RHS to comply with matrix regularization H = A + nu*B*B' */
1345: /* Add q = q + nu*B*b */
1346: if (jac->gkbnu) {
1347: nu = jac->gkbnu;
1348: VecScale(ilinkD->x,jac->gkbnu);
1349: MatMultAdd(jac->B,ilinkD->x,ilinkA->x,ilinkA->x); /* q = q + nu*B*b */
1350: } else {
1351: /* Situation when no augmented Lagrangian is used. Then we set inner */
1352: /* matrix N = I in [Ar13], and thus nu = 1. */
1353: nu = 1;
1354: }
1356: /* Transform rhs from [q,tilde{b}] to [0,b] */
1357: PetscLogEventBegin(ilinkA->event,ksp,ilinkA->x,ilinkA->y,NULL);
1358: KSPSolve(ksp,ilinkA->x,ilinkA->y);
1359: KSPCheckSolve(ksp,pc,ilinkA->y);
1360: PetscLogEventEnd(ilinkA->event,ksp,ilinkA->x,ilinkA->y,NULL);
1361: MatMultHermitianTranspose(jac->B,ilinkA->y,work1);
1362: VecAXPBY(work1,1.0/nu,-1.0,ilinkD->x); /* c = b - B'*x */
1364: /* First step of algorithm */
1365: VecNorm(work1,NORM_2,&beta); /* beta = sqrt(nu*c'*c)*/
1366: KSPCheckDot(ksp,beta);
1367: beta = PetscSqrtScalar(nu)*beta;
1368: VecAXPBY(v,nu/beta,0.0,work1); /* v = nu/beta *c */
1369: MatMult(jac->B,v,work2); /* u = H^{-1}*B*v */
1370: PetscLogEventBegin(ilinkA->event,ksp,work2,u,NULL);
1371: KSPSolve(ksp,work2,u);
1372: KSPCheckSolve(ksp,pc,u);
1373: PetscLogEventEnd(ilinkA->event,ksp,work2,u,NULL);
1374: MatMult(jac->H,u,Hu); /* alpha = u'*H*u */
1375: VecDot(Hu,u,&alpha);
1376: KSPCheckDot(ksp,alpha);
1377: if (PetscRealPart(alpha) <= 0.0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_NOT_CONVERGED,"GKB preconditioner diverged, H is not positive definite");
1378: alpha = PetscSqrtScalar(PetscAbsScalar(alpha));
1379: VecScale(u,1.0/alpha);
1380: VecAXPBY(d,1.0/alpha,0.0,v); /* v = nu/beta *c */
1382: z = beta/alpha;
1383: vecz[1] = z;
1385: /* Computation of first iterate x(1) and p(1) */
1386: VecAXPY(ilinkA->y,z,u);
1387: VecCopy(d,ilinkD->y);
1388: VecScale(ilinkD->y,-z);
1390: iterGKB = 1; lowbnd = 2*jac->gkbtol;
1391: if (jac->gkbmonitor) {
1392: PetscViewerASCIIPrintf(jac->gkbviewer,"%3D GKB Lower bound estimate %14.12e\n",iterGKB,lowbnd);
1393: }
1395: while (iterGKB < jac->gkbmaxit && lowbnd > jac->gkbtol) {
1396: iterGKB += 1;
1397: MatMultHermitianTranspose(jac->B,u,work1); /* v <- nu*(B'*u-alpha/nu*v) */
1398: VecAXPBY(v,nu,-alpha,work1);
1399: VecNorm(v,NORM_2,&beta); /* beta = sqrt(nu)*v'*v */
1400: beta = beta/PetscSqrtScalar(nu);
1401: VecScale(v,1.0/beta);
1402: MatMult(jac->B,v,work2); /* u <- H^{-1}*(B*v-beta*H*u) */
1403: MatMult(jac->H,u,Hu);
1404: VecAXPY(work2,-beta,Hu);
1405: PetscLogEventBegin(ilinkA->event,ksp,work2,u,NULL);
1406: KSPSolve(ksp,work2,u);
1407: KSPCheckSolve(ksp,pc,u);
1408: PetscLogEventEnd(ilinkA->event,ksp,work2,u,NULL);
1409: MatMult(jac->H,u,Hu); /* alpha = u'*H*u */
1410: VecDot(Hu,u,&alpha);
1411: KSPCheckDot(ksp,alpha);
1412: if (PetscRealPart(alpha) <= 0.0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_NOT_CONVERGED,"GKB preconditioner diverged, H is not positive definite");
1413: alpha = PetscSqrtScalar(PetscAbsScalar(alpha));
1414: VecScale(u,1.0/alpha);
1416: z = -beta/alpha*z; /* z <- beta/alpha*z */
1417: vecz[0] = z;
1419: /* Computation of new iterate x(i+1) and p(i+1) */
1420: VecAXPBY(d,1.0/alpha,-beta/alpha,v); /* d = (v-beta*d)/alpha */
1421: VecAXPY(ilinkA->y,z,u); /* r = r + z*u */
1422: VecAXPY(ilinkD->y,-z,d); /* p = p - z*d */
1423: MatMult(jac->H,ilinkA->y,Hu); /* ||u||_H = u'*H*u */
1424: VecDot(Hu,ilinkA->y,&nrmz2);
1426: /* Compute Lower Bound estimate */
1427: if (iterGKB > jac->gkbdelay) {
1428: lowbnd = 0.0;
1429: for (j=0; j<jac->gkbdelay; j++) {
1430: lowbnd += PetscAbsScalar(vecz[j]*vecz[j]);
1431: }
1432: lowbnd = PetscSqrtScalar(lowbnd/PetscAbsScalar(nrmz2));
1433: }
1435: for (j=0; j<jac->gkbdelay-1; j++) {
1436: vecz[jac->gkbdelay-j-1] = vecz[jac->gkbdelay-j-2];
1437: }
1438: if (jac->gkbmonitor) {
1439: PetscViewerASCIIPrintf(jac->gkbviewer,"%3D GKB Lower bound estimate %14.12e\n",iterGKB,lowbnd);
1440: }
1441: }
1443: VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
1444: VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);
1445: VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
1446: VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);
1448: return(0);
1449: }
1452: #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
1453: (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
1454: VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
1455: PetscLogEventBegin(ilink->event,ilink->ksp,ilink->y,ilink->x,NULL) || \
1456: KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
1457: KSPCheckSolve(ilink->ksp,pc,ilink->x) || \
1458: PetscLogEventEnd(ilink->event,ilink->ksp,ilink->y,ilink->x,NULL) || \
1459: VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
1460: VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
1462: static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
1463: {
1464: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1465: PetscErrorCode ierr;
1466: PC_FieldSplitLink ilink = jac->head;
1467: PetscInt bs;
1470: if (jac->type == PC_COMPOSITE_ADDITIVE) {
1471: if (jac->defaultsplit) {
1472: VecGetBlockSize(x,&bs);
1473: if (jac->bs > 0 && bs != jac->bs) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Blocksize of x vector %D does not match fieldsplit blocksize %D",bs,jac->bs);
1474: VecGetBlockSize(y,&bs);
1475: if (jac->bs > 0 && bs != jac->bs) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Blocksize of y vector %D does not match fieldsplit blocksize %D",bs,jac->bs);
1476: VecStrideGatherAll(x,jac->x,INSERT_VALUES);
1477: while (ilink) {
1478: PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1479: KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);
1480: KSPCheckSolve(ilink->ksp,pc,ilink->y);
1481: PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);
1482: ilink = ilink->next;
1483: }
1484: VecStrideScatterAll(jac->y,y,INSERT_VALUES);
1485: } else {
1486: VecSet(y,0.0);
1487: while (ilink) {
1488: FieldSplitSplitSolveAddTranspose(ilink,x,y);
1489: ilink = ilink->next;
1490: }
1491: }
1492: } else {
1493: if (!jac->w1) {
1494: VecDuplicate(x,&jac->w1);
1495: VecDuplicate(x,&jac->w2);
1496: }
1497: VecSet(y,0.0);
1498: if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
1499: FieldSplitSplitSolveAddTranspose(ilink,x,y);
1500: while (ilink->next) {
1501: ilink = ilink->next;
1502: MatMultTranspose(pc->mat,y,jac->w1);
1503: VecWAXPY(jac->w2,-1.0,jac->w1,x);
1504: FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);
1505: }
1506: while (ilink->previous) {
1507: ilink = ilink->previous;
1508: MatMultTranspose(pc->mat,y,jac->w1);
1509: VecWAXPY(jac->w2,-1.0,jac->w1,x);
1510: FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);
1511: }
1512: } else {
1513: while (ilink->next) { /* get to last entry in linked list */
1514: ilink = ilink->next;
1515: }
1516: FieldSplitSplitSolveAddTranspose(ilink,x,y);
1517: while (ilink->previous) {
1518: ilink = ilink->previous;
1519: MatMultTranspose(pc->mat,y,jac->w1);
1520: VecWAXPY(jac->w2,-1.0,jac->w1,x);
1521: FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);
1522: }
1523: }
1524: }
1525: return(0);
1526: }
1528: static PetscErrorCode PCReset_FieldSplit(PC pc)
1529: {
1530: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1531: PetscErrorCode ierr;
1532: PC_FieldSplitLink ilink = jac->head,next;
1535: while (ilink) {
1536: KSPDestroy(&ilink->ksp);
1537: VecDestroy(&ilink->x);
1538: VecDestroy(&ilink->y);
1539: VecDestroy(&ilink->z);
1540: VecScatterDestroy(&ilink->sctx);
1541: ISDestroy(&ilink->is);
1542: ISDestroy(&ilink->is_col);
1543: PetscFree(ilink->splitname);
1544: PetscFree(ilink->fields);
1545: PetscFree(ilink->fields_col);
1546: next = ilink->next;
1547: PetscFree(ilink);
1548: ilink = next;
1549: }
1550: jac->head = NULL;
1551: PetscFree2(jac->x,jac->y);
1552: if (jac->mat && jac->mat != jac->pmat) {
1553: MatDestroyMatrices(jac->nsplits,&jac->mat);
1554: } else if (jac->mat) {
1555: jac->mat = NULL;
1556: }
1557: if (jac->pmat) {MatDestroyMatrices(jac->nsplits,&jac->pmat);}
1558: if (jac->Afield) {MatDestroyMatrices(jac->nsplits,&jac->Afield);}
1559: jac->nsplits = 0;
1560: VecDestroy(&jac->w1);
1561: VecDestroy(&jac->w2);
1562: MatDestroy(&jac->schur);
1563: MatDestroy(&jac->schurp);
1564: MatDestroy(&jac->schur_user);
1565: KSPDestroy(&jac->kspschur);
1566: KSPDestroy(&jac->kspupper);
1567: MatDestroy(&jac->B);
1568: MatDestroy(&jac->C);
1569: MatDestroy(&jac->H);
1570: VecDestroy(&jac->u);
1571: VecDestroy(&jac->v);
1572: VecDestroy(&jac->Hu);
1573: VecDestroy(&jac->d);
1574: PetscFree(jac->vecz);
1575: PetscViewerDestroy(&jac->gkbviewer);
1576: jac->isrestrict = PETSC_FALSE;
1577: return(0);
1578: }
1580: static PetscErrorCode PCDestroy_FieldSplit(PC pc)
1581: {
1582: PetscErrorCode ierr;
1585: PCReset_FieldSplit(pc);
1586: PetscFree(pc->data);
1587: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSchurGetSubKSP_C",NULL);
1588: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",NULL);
1589: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",NULL);
1590: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",NULL);
1591: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",NULL);
1592: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",NULL);
1593: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",NULL);
1594: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",NULL);
1595: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",NULL);
1596: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitRestrictIS_C",NULL);
1597: return(0);
1598: }
1600: static PetscErrorCode PCSetFromOptions_FieldSplit(PetscOptionItems *PetscOptionsObject,PC pc)
1601: {
1602: PetscErrorCode ierr;
1603: PetscInt bs;
1604: PetscBool flg;
1605: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1606: PCCompositeType ctype;
1609: PetscOptionsHead(PetscOptionsObject,"FieldSplit options");
1610: PetscOptionsBool("-pc_fieldsplit_dm_splits","Whether to use DMCreateFieldDecomposition() for splits","PCFieldSplitSetDMSplits",jac->dm_splits,&jac->dm_splits,NULL);
1611: PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);
1612: if (flg) {
1613: PCFieldSplitSetBlockSize(pc,bs);
1614: }
1615: jac->diag_use_amat = pc->useAmat;
1616: PetscOptionsBool("-pc_fieldsplit_diag_use_amat","Use Amat (not Pmat) to extract diagonal fieldsplit blocks", "PCFieldSplitSetDiagUseAmat",jac->diag_use_amat,&jac->diag_use_amat,NULL);
1617: jac->offdiag_use_amat = pc->useAmat;
1618: PetscOptionsBool("-pc_fieldsplit_off_diag_use_amat","Use Amat (not Pmat) to extract off-diagonal fieldsplit blocks", "PCFieldSplitSetOffDiagUseAmat",jac->offdiag_use_amat,&jac->offdiag_use_amat,NULL);
1619: PetscOptionsBool("-pc_fieldsplit_detect_saddle_point","Form 2-way split by detecting zero diagonal entries", "PCFieldSplitSetDetectSaddlePoint",jac->detect,&jac->detect,NULL);
1620: PCFieldSplitSetDetectSaddlePoint(pc,jac->detect); /* Sets split type and Schur PC type */
1621: PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);
1622: if (flg) {
1623: PCFieldSplitSetType(pc,ctype);
1624: }
1625: /* Only setup fields once */
1626: if ((jac->bs > 0) && (jac->nsplits == 0)) {
1627: /* only allow user to set fields from command line if bs is already known.
1628: otherwise user can set them in PCFieldSplitSetDefaults() */
1629: PCFieldSplitSetRuntimeSplits_Private(pc);
1630: if (jac->splitdefined) {PetscInfo(pc,"Splits defined using the options database\n");}
1631: }
1632: if (jac->type == PC_COMPOSITE_SCHUR) {
1633: PetscOptionsGetEnum(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);
1634: if (flg) {PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");}
1635: PetscOptionsEnum("-pc_fieldsplit_schur_fact_type","Which off-diagonal parts of the block factorization to use","PCFieldSplitSetSchurFactType",PCFieldSplitSchurFactTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,NULL);
1636: PetscOptionsEnum("-pc_fieldsplit_schur_precondition","How to build preconditioner for Schur complement","PCFieldSplitSetSchurPre",PCFieldSplitSchurPreTypes,(PetscEnum)jac->schurpre,(PetscEnum*)&jac->schurpre,NULL);
1637: PetscOptionsScalar("-pc_fieldsplit_schur_scale","Scale Schur complement","PCFieldSplitSetSchurScale",jac->schurscale,&jac->schurscale,NULL);
1638: } else if (jac->type == PC_COMPOSITE_GKB) {
1639: PetscOptionsReal("-pc_fieldsplit_gkb_tol","The tolerance for the lower bound stopping criterion","PCFieldSplitGKBTol",jac->gkbtol,&jac->gkbtol,NULL);
1640: PetscOptionsInt("-pc_fieldsplit_gkb_delay","The delay value for lower bound criterion","PCFieldSplitGKBDelay",jac->gkbdelay,&jac->gkbdelay,NULL);
1641: PetscOptionsReal("-pc_fieldsplit_gkb_nu","Parameter in augmented Lagrangian approach","PCFieldSplitGKBNu",jac->gkbnu,&jac->gkbnu,NULL);
1642: if (jac->gkbnu < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nu cannot be less than 0: value %f",jac->gkbnu);
1643: PetscOptionsInt("-pc_fieldsplit_gkb_maxit","Maximum allowed number of iterations","PCFieldSplitGKBMaxit",jac->gkbmaxit,&jac->gkbmaxit,NULL);
1644: PetscOptionsBool("-pc_fieldsplit_gkb_monitor","Prints number of GKB iterations and error","PCFieldSplitGKB",jac->gkbmonitor,&jac->gkbmonitor,NULL);
1645: }
1646: PetscOptionsTail();
1647: return(0);
1648: }
1650: /*------------------------------------------------------------------------------------*/
1652: static PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
1653: {
1654: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1655: PetscErrorCode ierr;
1656: PC_FieldSplitLink ilink,next = jac->head;
1657: char prefix[128];
1658: PetscInt i;
1661: if (jac->splitdefined) {
1662: PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);
1663: return(0);
1664: }
1665: for (i=0; i<n; i++) {
1666: if (fields[i] >= jac->bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs);
1667: if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
1668: }
1669: PetscNew(&ilink);
1670: if (splitname) {
1671: PetscStrallocpy(splitname,&ilink->splitname);
1672: } else {
1673: PetscMalloc1(3,&ilink->splitname);
1674: PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);
1675: }
1676: ilink->event = jac->nsplits < 5 ? KSP_Solve_FS_0 + jac->nsplits : KSP_Solve_FS_0 + 4; /* Any split great than 4 gets logged in the 4th split */
1677: PetscMalloc1(n,&ilink->fields);
1678: PetscArraycpy(ilink->fields,fields,n);
1679: PetscMalloc1(n,&ilink->fields_col);
1680: PetscArraycpy(ilink->fields_col,fields_col,n);
1682: ilink->nfields = n;
1683: ilink->next = NULL;
1684: KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);
1685: KSPSetErrorIfNotConverged(ilink->ksp,pc->erroriffailure);
1686: PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);
1687: KSPSetType(ilink->ksp,KSPPREONLY);
1688: PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);
1690: PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);
1691: KSPSetOptionsPrefix(ilink->ksp,prefix);
1693: if (!next) {
1694: jac->head = ilink;
1695: ilink->previous = NULL;
1696: } else {
1697: while (next->next) {
1698: next = next->next;
1699: }
1700: next->next = ilink;
1701: ilink->previous = next;
1702: }
1703: jac->nsplits++;
1704: return(0);
1705: }
1707: static PetscErrorCode PCFieldSplitSchurGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
1708: {
1709: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1713: *subksp = NULL;
1714: if (n) *n = 0;
1715: if (jac->type == PC_COMPOSITE_SCHUR) {
1716: PetscInt nn;
1718: if (!jac->schur) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call KSPSetUp() or PCSetUp() before calling PCFieldSplitSchurGetSubKSP()");
1719: if (jac->nsplits != 2) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Unexpected number of splits %D != 2",jac->nsplits);
1720: nn = jac->nsplits + (jac->kspupper != jac->head->ksp ? 1 : 0);
1721: PetscMalloc1(nn,subksp);
1722: (*subksp)[0] = jac->head->ksp;
1723: (*subksp)[1] = jac->kspschur;
1724: if (jac->kspupper != jac->head->ksp) (*subksp)[2] = jac->kspupper;
1725: if (n) *n = nn;
1726: }
1727: return(0);
1728: }
1730: static PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
1731: {
1732: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1736: if (!jac->schur) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call KSPSetUp() or PCSetUp() before calling PCFieldSplitGetSubKSP()");
1737: PetscMalloc1(jac->nsplits,subksp);
1738: MatSchurComplementGetKSP(jac->schur,*subksp);
1740: (*subksp)[1] = jac->kspschur;
1741: if (n) *n = jac->nsplits;
1742: return(0);
1743: }
1745: static PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
1746: {
1747: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1748: PetscErrorCode ierr;
1749: PetscInt cnt = 0;
1750: PC_FieldSplitLink ilink = jac->head;
1753: PetscMalloc1(jac->nsplits,subksp);
1754: while (ilink) {
1755: (*subksp)[cnt++] = ilink->ksp;
1756: ilink = ilink->next;
1757: }
1758: if (cnt != jac->nsplits) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number of splits in linked list %D does not match number in object %D",cnt,jac->nsplits);
1759: if (n) *n = jac->nsplits;
1760: return(0);
1761: }
1763: /*@C
1764: PCFieldSplitRestrictIS - Restricts the fieldsplit ISs to be within a given IS.
1766: Input Parameters:
1767: + pc - the preconditioner context
1768: - is - the index set that defines the indices to which the fieldsplit is to be restricted
1770: Level: advanced
1772: @*/
1773: PetscErrorCode PCFieldSplitRestrictIS(PC pc,IS isy)
1774: {
1780: PetscTryMethod(pc,"PCFieldSplitRestrictIS_C",(PC,IS),(pc,isy));
1781: return(0);
1782: }
1785: static PetscErrorCode PCFieldSplitRestrictIS_FieldSplit(PC pc, IS isy)
1786: {
1787: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1788: PetscErrorCode ierr;
1789: PC_FieldSplitLink ilink = jac->head, next;
1790: PetscInt localsize,size,sizez,i;
1791: const PetscInt *ind, *indz;
1792: PetscInt *indc, *indcz;
1793: PetscBool flg;
1796: ISGetLocalSize(isy,&localsize);
1797: MPI_Scan(&localsize,&size,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)isy));
1798: size -= localsize;
1799: while(ilink) {
1800: IS isrl,isr;
1801: PC subpc;
1802: ISEmbed(ilink->is, isy, PETSC_TRUE, &isrl);
1803: ISGetLocalSize(isrl,&localsize);
1804: PetscMalloc1(localsize,&indc);
1805: ISGetIndices(isrl,&ind);
1806: PetscArraycpy(indc,ind,localsize);
1807: ISRestoreIndices(isrl,&ind);
1808: ISDestroy(&isrl);
1809: for (i=0; i<localsize; i++) *(indc+i) += size;
1810: ISCreateGeneral(PetscObjectComm((PetscObject)isy),localsize,indc,PETSC_OWN_POINTER,&isr);
1811: PetscObjectReference((PetscObject)isr);
1812: ISDestroy(&ilink->is);
1813: ilink->is = isr;
1814: PetscObjectReference((PetscObject)isr);
1815: ISDestroy(&ilink->is_col);
1816: ilink->is_col = isr;
1817: ISDestroy(&isr);
1818: KSPGetPC(ilink->ksp, &subpc);
1819: PetscObjectTypeCompare((PetscObject)subpc,PCFIELDSPLIT,&flg);
1820: if(flg) {
1821: IS iszl,isz;
1822: MPI_Comm comm;
1823: ISGetLocalSize(ilink->is,&localsize);
1824: comm = PetscObjectComm((PetscObject)ilink->is);
1825: ISEmbed(isy, ilink->is, PETSC_TRUE, &iszl);
1826: MPI_Scan(&localsize,&sizez,1,MPIU_INT,MPI_SUM,comm);
1827: sizez -= localsize;
1828: ISGetLocalSize(iszl,&localsize);
1829: PetscMalloc1(localsize,&indcz);
1830: ISGetIndices(iszl,&indz);
1831: PetscArraycpy(indcz,indz,localsize);
1832: ISRestoreIndices(iszl,&indz);
1833: ISDestroy(&iszl);
1834: for (i=0; i<localsize; i++) *(indcz+i) += sizez;
1835: ISCreateGeneral(comm,localsize,indcz,PETSC_OWN_POINTER,&isz);
1836: PCFieldSplitRestrictIS(subpc,isz);
1837: ISDestroy(&isz);
1838: }
1839: next = ilink->next;
1840: ilink = next;
1841: }
1842: jac->isrestrict = PETSC_TRUE;
1843: return(0);
1844: }
1846: static PetscErrorCode PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
1847: {
1848: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1849: PetscErrorCode ierr;
1850: PC_FieldSplitLink ilink, next = jac->head;
1851: char prefix[128];
1854: if (jac->splitdefined) {
1855: PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);
1856: return(0);
1857: }
1858: PetscNew(&ilink);
1859: if (splitname) {
1860: PetscStrallocpy(splitname,&ilink->splitname);
1861: } else {
1862: PetscMalloc1(8,&ilink->splitname);
1863: PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);
1864: }
1865: ilink->event = jac->nsplits < 5 ? KSP_Solve_FS_0 + jac->nsplits : KSP_Solve_FS_0 + 4; /* Any split great than 4 gets logged in the 4th split */
1866: PetscObjectReference((PetscObject)is);
1867: ISDestroy(&ilink->is);
1868: ilink->is = is;
1869: PetscObjectReference((PetscObject)is);
1870: ISDestroy(&ilink->is_col);
1871: ilink->is_col = is;
1872: ilink->next = NULL;
1873: KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);
1874: KSPSetErrorIfNotConverged(ilink->ksp,pc->erroriffailure);
1875: PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);
1876: KSPSetType(ilink->ksp,KSPPREONLY);
1877: PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);
1879: PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);
1880: KSPSetOptionsPrefix(ilink->ksp,prefix);
1882: if (!next) {
1883: jac->head = ilink;
1884: ilink->previous = NULL;
1885: } else {
1886: while (next->next) {
1887: next = next->next;
1888: }
1889: next->next = ilink;
1890: ilink->previous = next;
1891: }
1892: jac->nsplits++;
1893: return(0);
1894: }
1896: /*@C
1897: PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
1899: Logically Collective on PC
1901: Input Parameters:
1902: + pc - the preconditioner context
1903: . splitname - name of this split, if NULL the number of the split is used
1904: . n - the number of fields in this split
1905: - fields - the fields in this split
1907: Level: intermediate
1909: Notes:
1910: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
1912: The PCFieldSplitSetFields() is for defining fields as strided blocks. For example, if the block
1913: size is three then one can define a field as 0, or 1 or 2 or 0,1 or 0,2 or 1,2 which mean
1914: 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1915: where the numbered entries indicate what is in the field.
1917: This function is called once per split (it creates a new split each time). Solve options
1918: for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1920: Developer Note: This routine does not actually create the IS representing the split, that is delayed
1921: until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be
1922: available when this routine is called.
1924: .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
1926: @*/
1927: PetscErrorCode PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
1928: {
1934: if (n < 1) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1936: PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt*,const PetscInt*),(pc,splitname,n,fields,fields_col));
1937: return(0);
1938: }
1940: /*@
1941: PCFieldSplitSetDiagUseAmat - set flag indicating whether to extract diagonal blocks from Amat (rather than Pmat)
1943: Logically Collective on PC
1945: Input Parameters:
1946: + pc - the preconditioner object
1947: - flg - boolean flag indicating whether or not to use Amat to extract the diagonal blocks from
1949: Options Database:
1950: . -pc_fieldsplit_diag_use_amat
1952: Level: intermediate
1954: .seealso: PCFieldSplitGetDiagUseAmat(), PCFieldSplitSetOffDiagUseAmat(), PCFIELDSPLIT
1956: @*/
1957: PetscErrorCode PCFieldSplitSetDiagUseAmat(PC pc,PetscBool flg)
1958: {
1959: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1960: PetscBool isfs;
1965: PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
1966: if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1967: jac->diag_use_amat = flg;
1968: return(0);
1969: }
1971: /*@
1972: PCFieldSplitGetDiagUseAmat - get the flag indicating whether to extract diagonal blocks from Amat (rather than Pmat)
1974: Logically Collective on PC
1976: Input Parameters:
1977: . pc - the preconditioner object
1979: Output Parameters:
1980: . flg - boolean flag indicating whether or not to use Amat to extract the diagonal blocks from
1983: Level: intermediate
1985: .seealso: PCFieldSplitSetDiagUseAmat(), PCFieldSplitGetOffDiagUseAmat(), PCFIELDSPLIT
1987: @*/
1988: PetscErrorCode PCFieldSplitGetDiagUseAmat(PC pc,PetscBool *flg)
1989: {
1990: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1991: PetscBool isfs;
1997: PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
1998: if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1999: *flg = jac->diag_use_amat;
2000: return(0);
2001: }
2003: /*@
2004: PCFieldSplitSetOffDiagUseAmat - set flag indicating whether to extract off-diagonal blocks from Amat (rather than Pmat)
2006: Logically Collective on PC
2008: Input Parameters:
2009: + pc - the preconditioner object
2010: - flg - boolean flag indicating whether or not to use Amat to extract the off-diagonal blocks from
2012: Options Database:
2013: . -pc_fieldsplit_off_diag_use_amat
2015: Level: intermediate
2017: .seealso: PCFieldSplitGetOffDiagUseAmat(), PCFieldSplitSetDiagUseAmat(), PCFIELDSPLIT
2019: @*/
2020: PetscErrorCode PCFieldSplitSetOffDiagUseAmat(PC pc,PetscBool flg)
2021: {
2022: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2023: PetscBool isfs;
2028: PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
2029: if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
2030: jac->offdiag_use_amat = flg;
2031: return(0);
2032: }
2034: /*@
2035: PCFieldSplitGetOffDiagUseAmat - get the flag indicating whether to extract off-diagonal blocks from Amat (rather than Pmat)
2037: Logically Collective on PC
2039: Input Parameters:
2040: . pc - the preconditioner object
2042: Output Parameters:
2043: . flg - boolean flag indicating whether or not to use Amat to extract the off-diagonal blocks from
2046: Level: intermediate
2048: .seealso: PCFieldSplitSetOffDiagUseAmat(), PCFieldSplitGetDiagUseAmat(), PCFIELDSPLIT
2050: @*/
2051: PetscErrorCode PCFieldSplitGetOffDiagUseAmat(PC pc,PetscBool *flg)
2052: {
2053: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2054: PetscBool isfs;
2060: PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
2061: if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
2062: *flg = jac->offdiag_use_amat;
2063: return(0);
2064: }
2068: /*@C
2069: PCFieldSplitSetIS - Sets the exact elements for field
2071: Logically Collective on PC
2073: Input Parameters:
2074: + pc - the preconditioner context
2075: . splitname - name of this split, if NULL the number of the split is used
2076: - is - the index set that defines the vector elements in this field
2079: Notes:
2080: Use PCFieldSplitSetFields(), for fields defined by strided types.
2082: This function is called once per split (it creates a new split each time). Solve options
2083: for this split will be available under the prefix -fieldsplit_SPLITNAME_.
2085: Level: intermediate
2087: .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
2089: @*/
2090: PetscErrorCode PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
2091: {
2098: PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));
2099: return(0);
2100: }
2102: /*@C
2103: PCFieldSplitGetIS - Retrieves the elements for a field as an IS
2105: Logically Collective on PC
2107: Input Parameters:
2108: + pc - the preconditioner context
2109: - splitname - name of this split
2111: Output Parameter:
2112: - is - the index set that defines the vector elements in this field, or NULL if the field is not found
2114: Level: intermediate
2116: .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS()
2118: @*/
2119: PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
2120: {
2127: {
2128: PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;
2129: PC_FieldSplitLink ilink = jac->head;
2130: PetscBool found;
2132: *is = NULL;
2133: while (ilink) {
2134: PetscStrcmp(ilink->splitname, splitname, &found);
2135: if (found) {
2136: *is = ilink->is;
2137: break;
2138: }
2139: ilink = ilink->next;
2140: }
2141: }
2142: return(0);
2143: }
2145: /*@C
2146: PCFieldSplitGetISByIndex - Retrieves the elements for a given index field as an IS
2148: Logically Collective on PC
2150: Input Parameters:
2151: + pc - the preconditioner context
2152: - index - index of this split
2154: Output Parameter:
2155: - is - the index set that defines the vector elements in this field
2157: Level: intermediate
2159: .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitGetIS(), PCFieldSplitSetIS()
2161: @*/
2162: PetscErrorCode PCFieldSplitGetISByIndex(PC pc,PetscInt index,IS *is)
2163: {
2167: if (index < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",index);
2170: {
2171: PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;
2172: PC_FieldSplitLink ilink = jac->head;
2173: PetscInt i = 0;
2174: if (index >= jac->nsplits) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",index,jac->nsplits);
2176: while (i < index) {
2177: ilink = ilink->next;
2178: ++i;
2179: }
2180: PCFieldSplitGetIS(pc,ilink->splitname,is);
2181: }
2182: return(0);
2183: }
2185: /*@
2186: PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
2187: fieldsplit preconditioner. If not set the matrix block size is used.
2189: Logically Collective on PC
2191: Input Parameters:
2192: + pc - the preconditioner context
2193: - bs - the block size
2195: Level: intermediate
2197: .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
2199: @*/
2200: PetscErrorCode PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
2201: {
2207: PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));
2208: return(0);
2209: }
2211: /*@C
2212: PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
2214: Collective on KSP
2216: Input Parameter:
2217: . pc - the preconditioner context
2219: Output Parameters:
2220: + n - the number of splits
2221: - subksp - the array of KSP contexts
2223: Note:
2224: After PCFieldSplitGetSubKSP() the array of KSPs is to be freed by the user with PetscFree()
2225: (not the KSP just the array that contains them).
2227: You must call PCSetUp() before calling PCFieldSplitGetSubKSP().
2229: If the fieldsplit is of type PC_COMPOSITE_SCHUR, it returns the KSP object used inside the
2230: Schur complement and the KSP object used to iterate over the Schur complement.
2231: To access all the KSP objects used in PC_COMPOSITE_SCHUR, use PCFieldSplitSchurGetSubKSP().
2233: If the fieldsplit is of type PC_COMPOSITE_GKB, it returns the KSP object used to solve the
2234: inner linear system defined by the matrix H in each loop.
2236: Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs.
2237: You can call PCFieldSplitGetSubKSP(pc,n,PETSC_NULL_KSP,ierr) to determine how large the
2238: KSP array must be.
2241: Level: advanced
2243: .seealso: PCFIELDSPLIT
2244: @*/
2245: PetscErrorCode PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
2246: {
2252: PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));
2253: return(0);
2254: }
2256: /*@C
2257: PCFieldSplitSchurGetSubKSP - Gets the KSP contexts used inside the Schur complement based PCFIELDSPLIT
2259: Collective on KSP
2261: Input Parameter:
2262: . pc - the preconditioner context
2264: Output Parameters:
2265: + n - the number of splits
2266: - subksp - the array of KSP contexts
2268: Note:
2269: After PCFieldSplitSchurGetSubKSP() the array of KSPs is to be freed by the user with PetscFree()
2270: (not the KSP just the array that contains them).
2272: You must call PCSetUp() before calling PCFieldSplitSchurGetSubKSP().
2274: If the fieldsplit type is of type PC_COMPOSITE_SCHUR, it returns (in order)
2275: - the KSP used for the (1,1) block
2276: - the KSP used for the Schur complement (not the one used for the interior Schur solver)
2277: - the KSP used for the (1,1) block in the upper triangular factor (if different from that of the (1,1) block).
2279: It returns a null array if the fieldsplit is not of type PC_COMPOSITE_SCHUR; in this case, you should use PCFieldSplitGetSubKSP().
2281: Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs.
2282: You can call PCFieldSplitSchurGetSubKSP(pc,n,PETSC_NULL_KSP,ierr) to determine how large the
2283: KSP array must be.
2285: Level: advanced
2287: .seealso: PCFIELDSPLIT
2288: @*/
2289: PetscErrorCode PCFieldSplitSchurGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
2290: {
2296: PetscUseMethod(pc,"PCFieldSplitSchurGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));
2297: return(0);
2298: }
2300: /*@
2301: PCFieldSplitSetSchurPre - Indicates what operator is used to construct the preconditioner for the Schur complement.
2302: A11 matrix. Otherwise no preconditioner is used.
2304: Collective on PC
2306: Input Parameters:
2307: + pc - the preconditioner context
2308: . ptype - which matrix to use for preconditioning the Schur complement: PC_FIELDSPLIT_SCHUR_PRE_A11 (default), PC_FIELDSPLIT_SCHUR_PRE_SELF, PC_FIELDSPLIT_SCHUR_PRE_USER
2309: PC_FIELDSPLIT_SCHUR_PRE_SELFP, and PC_FIELDSPLIT_SCHUR_PRE_FULL
2310: - userpre - matrix to use for preconditioning, or NULL
2312: Options Database:
2313: . -pc_fieldsplit_schur_precondition <self,selfp,user,a11,full> - default is a11. See notes for meaning of various arguments
2315: Notes:
2316: $ If ptype is
2317: $ a11 then the preconditioner for the Schur complement is generated from the block diagonal part of the preconditioner
2318: $ matrix associated with the Schur complement (i.e. A11), not the Schur complement matrix
2319: $ self the preconditioner for the Schur complement is generated from the symbolic representation of the Schur complement matrix:
2320: $ The only preconditioner that currently works with this symbolic respresentation matrix object is the PCLSC
2321: $ preconditioner
2322: $ user then the preconditioner for the Schur complement is generated from the user provided matrix (pre argument
2323: $ to this function).
2324: $ selfp then the preconditioning for the Schur complement is generated from an explicitly-assembled approximation Sp = A11 - A10 inv(diag(A00)) A01
2325: $ This is only a good preconditioner when diag(A00) is a good preconditioner for A00. Optionally, A00 can be
2326: $ lumped before extracting the diagonal using the additional option -fieldsplit_1_mat_schur_complement_ainv_type lump
2327: $ full then the preconditioner for the Schur complement is generated from the exact Schur complement matrix representation computed internally by PCFIELDSPLIT (this is expensive)
2328: $ useful mostly as a test that the Schur complement approach can work for your problem
2330: When solving a saddle point problem, where the A11 block is identically zero, using a11 as the ptype only makes sense
2331: with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
2332: -fieldsplit_1_pc_type lsc which uses the least squares commutator to compute a preconditioner for the Schur complement.
2334: Level: intermediate
2336: .seealso: PCFieldSplitGetSchurPre(), PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType,
2337: MatSchurComplementSetAinvType(), PCLSC
2339: @*/
2340: PetscErrorCode PCFieldSplitSetSchurPre(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
2341: {
2346: PetscTryMethod(pc,"PCFieldSplitSetSchurPre_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));
2347: return(0);
2348: }
2350: PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) {return PCFieldSplitSetSchurPre(pc,ptype,pre);} /* Deprecated name */
2352: /*@
2353: PCFieldSplitGetSchurPre - For Schur complement fieldsplit, determine how the Schur complement will be
2354: preconditioned. See PCFieldSplitSetSchurPre() for details.
2356: Logically Collective on PC
2358: Input Parameters:
2359: . pc - the preconditioner context
2361: Output Parameters:
2362: + ptype - which matrix to use for preconditioning the Schur complement: PC_FIELDSPLIT_SCHUR_PRE_A11, PC_FIELDSPLIT_SCHUR_PRE_SELF, PC_FIELDSPLIT_PRE_USER
2363: - userpre - matrix to use for preconditioning (with PC_FIELDSPLIT_PRE_USER), or NULL
2365: Level: intermediate
2367: .seealso: PCFieldSplitSetSchurPre(), PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC
2369: @*/
2370: PetscErrorCode PCFieldSplitGetSchurPre(PC pc,PCFieldSplitSchurPreType *ptype,Mat *pre)
2371: {
2376: PetscUseMethod(pc,"PCFieldSplitGetSchurPre_C",(PC,PCFieldSplitSchurPreType*,Mat*),(pc,ptype,pre));
2377: return(0);
2378: }
2380: /*@
2381: PCFieldSplitSchurGetS - extract the MatSchurComplement object used by this PC in case it needs to be configured separately
2383: Not collective
2385: Input Parameter:
2386: . pc - the preconditioner context
2388: Output Parameter:
2389: . S - the Schur complement matrix
2391: Notes:
2392: This matrix should not be destroyed using MatDestroy(); rather, use PCFieldSplitSchurRestoreS().
2394: Level: advanced
2396: .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSchurPreType, PCFieldSplitSetSchurPre(), MatSchurComplement, PCFieldSplitSchurRestoreS()
2398: @*/
2399: PetscErrorCode PCFieldSplitSchurGetS(PC pc,Mat *S)
2400: {
2402: const char* t;
2403: PetscBool isfs;
2404: PC_FieldSplit *jac;
2408: PetscObjectGetType((PetscObject)pc,&t);
2409: PetscStrcmp(t,PCFIELDSPLIT,&isfs);
2410: if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PC of type PCFIELDSPLIT, got %s instead",t);
2411: jac = (PC_FieldSplit*)pc->data;
2412: if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PCFIELDSPLIT of type SCHUR, got %D instead",jac->type);
2413: if (S) *S = jac->schur;
2414: return(0);
2415: }
2417: /*@
2418: PCFieldSplitSchurRestoreS - restores the MatSchurComplement object used by this PC
2420: Not collective
2422: Input Parameters:
2423: + pc - the preconditioner context
2424: - S - the Schur complement matrix
2426: Level: advanced
2428: .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSchurPreType, PCFieldSplitSetSchurPre(), MatSchurComplement, PCFieldSplitSchurGetS()
2430: @*/
2431: PetscErrorCode PCFieldSplitSchurRestoreS(PC pc,Mat *S)
2432: {
2434: const char* t;
2435: PetscBool isfs;
2436: PC_FieldSplit *jac;
2440: PetscObjectGetType((PetscObject)pc,&t);
2441: PetscStrcmp(t,PCFIELDSPLIT,&isfs);
2442: if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PC of type PCFIELDSPLIT, got %s instead",t);
2443: jac = (PC_FieldSplit*)pc->data;
2444: if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PCFIELDSPLIT of type SCHUR, got %D instead",jac->type);
2445: if (!S || *S != jac->schur) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"MatSchurComplement restored is not the same as gotten");
2446: return(0);
2447: }
2450: static PetscErrorCode PCFieldSplitSetSchurPre_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
2451: {
2452: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2456: jac->schurpre = ptype;
2457: if (ptype == PC_FIELDSPLIT_SCHUR_PRE_USER && pre) {
2458: MatDestroy(&jac->schur_user);
2459: jac->schur_user = pre;
2460: PetscObjectReference((PetscObject)jac->schur_user);
2461: }
2462: return(0);
2463: }
2465: static PetscErrorCode PCFieldSplitGetSchurPre_FieldSplit(PC pc,PCFieldSplitSchurPreType *ptype,Mat *pre)
2466: {
2467: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2470: *ptype = jac->schurpre;
2471: *pre = jac->schur_user;
2472: return(0);
2473: }
2475: /*@
2476: PCFieldSplitSetSchurFactType - sets which blocks of the approximate block factorization to retain in the preconditioner
2478: Collective on PC
2480: Input Parameters:
2481: + pc - the preconditioner context
2482: - ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default
2484: Options Database:
2485: . -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full
2488: Level: intermediate
2490: Notes:
2491: The FULL factorization is
2493: $ (A B) = (1 0) (A 0) (1 Ainv*B) = L D U
2494: $ (C E) (C*Ainv 1) (0 S) (0 1 )
2496: where S = E - C*Ainv*B. In practice, the full factorization is applied via block triangular solves with the grouping L*(D*U). UPPER uses D*U, LOWER uses L*D,
2497: and DIAG is the diagonal part with the sign of S flipped (because this makes the preconditioner positive definite for many formulations, thus allowing the use of KSPMINRES). Sign flipping of S can be turned off with PCFieldSplitSetSchurScale().
2499: $ If A and S are solved exactly
2500: $ *) FULL factorization is a direct solver.
2501: $ *) The preconditioned operator with LOWER or UPPER has all eigenvalues equal to 1 and minimal polynomial of degree 2, so KSPGMRES converges in 2 iterations.
2502: $ *) With DIAG, the preconditioned operator has three distinct nonzero eigenvalues and minimal polynomial of degree at most 4, so KSPGMRES converges in at most 4 iterations.
2504: If the iteration count is very low, consider using KSPFGMRES or KSPGCR which can use one less preconditioner
2505: application in this case. Note that the preconditioned operator may be highly non-normal, so such fast convergence may not be observed in practice.
2507: For symmetric problems in which A is positive definite and S is negative definite, DIAG can be used with KSPMINRES.
2509: Note that a flexible method like KSPFGMRES or KSPGCR must be used if the fieldsplit preconditioner is nonlinear (e.g. a few iterations of a Krylov method is used to solve with A or S).
2511: References:
2512: + 1. - Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000).
2513: - 2. - Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001).
2515: .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCFieldSplitSetSchurScale()
2516: @*/
2517: PetscErrorCode PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype)
2518: {
2523: PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));
2524: return(0);
2525: }
2527: static PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype)
2528: {
2529: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2532: jac->schurfactorization = ftype;
2533: return(0);
2534: }
2536: /*@
2537: PCFieldSplitSetSchurScale - Controls the sign flip of S for PC_FIELDSPLIT_SCHUR_FACT_DIAG.
2539: Collective on PC
2541: Input Parameters:
2542: + pc - the preconditioner context
2543: - scale - scaling factor for the Schur complement
2545: Options Database:
2546: . -pc_fieldsplit_schur_scale - default is -1.0
2548: Level: intermediate
2550: .seealso: PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurFactType, PCFieldSplitSetSchurScale()
2551: @*/
2552: PetscErrorCode PCFieldSplitSetSchurScale(PC pc,PetscScalar scale)
2553: {
2559: PetscTryMethod(pc,"PCFieldSplitSetSchurScale_C",(PC,PetscScalar),(pc,scale));
2560: return(0);
2561: }
2563: static PetscErrorCode PCFieldSplitSetSchurScale_FieldSplit(PC pc,PetscScalar scale)
2564: {
2565: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2568: jac->schurscale = scale;
2569: return(0);
2570: }
2572: /*@C
2573: PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
2575: Collective on KSP
2577: Input Parameter:
2578: . pc - the preconditioner context
2580: Output Parameters:
2581: + A00 - the (0,0) block
2582: . A01 - the (0,1) block
2583: . A10 - the (1,0) block
2584: - A11 - the (1,1) block
2586: Level: advanced
2588: .seealso: PCFIELDSPLIT
2589: @*/
2590: PetscErrorCode PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
2591: {
2592: PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;
2596: if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
2597: if (A00) *A00 = jac->pmat[0];
2598: if (A01) *A01 = jac->B;
2599: if (A10) *A10 = jac->C;
2600: if (A11) *A11 = jac->pmat[1];
2601: return(0);
2602: }
2604: /*@
2605: PCFieldSplitSetGKBTol - Sets the solver tolerance for the generalized Golub-Kahan bidiagonalization preconditioner.
2607: Collective on PC
2609: Notes:
2610: The generalized GKB algorithm uses a lower bound estimate of the error in energy norm as stopping criterion.
2611: It stops once the lower bound estimate undershoots the required solver tolerance. Although the actual error might be bigger than
2612: this estimate, the stopping criterion is satisfactory in practical cases [A13].
2614: [Ar13] Generalized Golub-Kahan bidiagonalization and stopping criteria, SIAM J. Matrix Anal. Appl., Vol. 34, No. 2, pp. 571-592, 2013.
2616: Input Parameters:
2617: + pc - the preconditioner context
2618: - tolerance - the solver tolerance
2620: Options Database:
2621: . -pc_fieldsplit_gkb_tol - default is 1e-5
2623: Level: intermediate
2625: .seealso: PCFIELDSPLIT, PCFieldSplitSetGKBDelay(), PCFieldSplitSetGKBNu(), PCFieldSplitSetGKBMaxit()
2626: @*/
2627: PetscErrorCode PCFieldSplitSetGKBTol(PC pc,PetscReal tolerance)
2628: {
2634: PetscTryMethod(pc,"PCFieldSplitSetGKBTol_C",(PC,PetscReal),(pc,tolerance));
2635: return(0);
2636: }
2638: static PetscErrorCode PCFieldSplitSetGKBTol_FieldSplit(PC pc,PetscReal tolerance)
2639: {
2640: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2643: jac->gkbtol = tolerance;
2644: return(0);
2645: }
2648: /*@
2649: PCFieldSplitSetGKBMaxit - Sets the maximum number of iterations for the generalized Golub-Kahan bidiagonalization
2650: preconditioner.
2652: Collective on PC
2654: Input Parameters:
2655: + pc - the preconditioner context
2656: - maxit - the maximum number of iterations
2658: Options Database:
2659: . -pc_fieldsplit_gkb_maxit - default is 100
2661: Level: intermediate
2663: .seealso: PCFIELDSPLIT, PCFieldSplitSetGKBDelay(), PCFieldSplitSetGKBTol(), PCFieldSplitSetGKBNu()
2664: @*/
2665: PetscErrorCode PCFieldSplitSetGKBMaxit(PC pc,PetscInt maxit)
2666: {
2672: PetscTryMethod(pc,"PCFieldSplitSetGKBMaxit_C",(PC,PetscInt),(pc,maxit));
2673: return(0);
2674: }
2676: static PetscErrorCode PCFieldSplitSetGKBMaxit_FieldSplit(PC pc,PetscInt maxit)
2677: {
2678: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2681: jac->gkbmaxit = maxit;
2682: return(0);
2683: }
2685: /*@
2686: PCFieldSplitSetGKBDelay - Sets the delay in the lower bound error estimate in the generalized Golub-Kahan bidiagonalization
2687: preconditioner.
2689: Collective on PC
2691: Notes:
2692: The algorithm uses a lower bound estimate of the error in energy norm as stopping criterion. The lower bound of the error ||u-u^k||_H
2693: is expressed as a truncated sum. The error at iteration k can only be measured at iteration (k + delay), and thus the algorithm needs
2694: at least (delay + 1) iterations to stop. For more details on the generalized Golub-Kahan bidiagonalization method and its lower bound stopping criterion, please refer to
2696: [Ar13] Generalized Golub-Kahan bidiagonalization and stopping criteria, SIAM J. Matrix Anal. Appl., Vol. 34, No. 2, pp. 571-592, 2013.
2698: Input Parameters:
2699: + pc - the preconditioner context
2700: - delay - the delay window in the lower bound estimate
2702: Options Database:
2703: . -pc_fieldsplit_gkb_delay - default is 5
2705: Level: intermediate
2707: .seealso: PCFIELDSPLIT, PCFieldSplitSetGKBNu(), PCFieldSplitSetGKBTol(), PCFieldSplitSetGKBMaxit()
2708: @*/
2709: PetscErrorCode PCFieldSplitSetGKBDelay(PC pc,PetscInt delay)
2710: {
2716: PetscTryMethod(pc,"PCFieldSplitSetGKBDelay_C",(PC,PetscInt),(pc,delay));
2717: return(0);
2718: }
2720: static PetscErrorCode PCFieldSplitSetGKBDelay_FieldSplit(PC pc,PetscInt delay)
2721: {
2722: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2725: jac->gkbdelay = delay;
2726: return(0);
2727: }
2729: /*@
2730: PCFieldSplitSetGKBNu - Sets the scalar value nu >= 0 in the transformation H = A00 + nu*A01*A01' of the (1,1) block in the Golub-Kahan bidiagonalization preconditioner.
2732: Collective on PC
2734: Notes:
2735: This shift is in general done to obtain better convergence properties for the outer loop of the algorithm. This is often achieved by chosing nu sufficiently big. However,
2736: if nu is chosen too big, the matrix H might be badly conditioned and the solution of the linear system Hx = b in the inner loop gets difficult. It is therefore
2737: necessary to find a good balance in between the convergence of the inner and outer loop.
2739: For nu = 0, no shift is done. In this case A00 has to be positive definite. The matrix N in [Ar13] is then chosen as identity.
2741: [Ar13] Generalized Golub-Kahan bidiagonalization and stopping criteria, SIAM J. Matrix Anal. Appl., Vol. 34, No. 2, pp. 571-592, 2013.
2743: Input Parameters:
2744: + pc - the preconditioner context
2745: - nu - the shift parameter
2747: Options Database:
2748: . -pc_fieldsplit_gkb_nu - default is 1
2750: Level: intermediate
2752: .seealso: PCFIELDSPLIT, PCFieldSplitSetGKBDelay(), PCFieldSplitSetGKBTol(), PCFieldSplitSetGKBMaxit()
2753: @*/
2754: PetscErrorCode PCFieldSplitSetGKBNu(PC pc,PetscReal nu)
2755: {
2761: PetscTryMethod(pc,"PCFieldSplitSetGKBNu_C",(PC,PetscReal),(pc,nu));
2762: return(0);
2763: }
2765: static PetscErrorCode PCFieldSplitSetGKBNu_FieldSplit(PC pc,PetscReal nu)
2766: {
2767: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2770: jac->gkbnu = nu;
2771: return(0);
2772: }
2775: static PetscErrorCode PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
2776: {
2777: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2781: jac->type = type;
2783: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",0);
2784: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",0);
2785: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",0);
2786: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",0);
2787: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurScale_C",0);
2788: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetGKBTol_C",0);
2789: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetGKBMaxit_C",0);
2790: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetGKBNu_C",0);
2791: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetGKBDelay_C",0);
2793: if (type == PC_COMPOSITE_SCHUR) {
2794: pc->ops->apply = PCApply_FieldSplit_Schur;
2795: pc->ops->view = PCView_FieldSplit_Schur;
2797: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit_Schur);
2798: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",PCFieldSplitSetSchurPre_FieldSplit);
2799: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",PCFieldSplitGetSchurPre_FieldSplit);
2800: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",PCFieldSplitSetSchurFactType_FieldSplit);
2801: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurScale_C",PCFieldSplitSetSchurScale_FieldSplit);
2802: } else if (type == PC_COMPOSITE_GKB){
2803: pc->ops->apply = PCApply_FieldSplit_GKB;
2804: pc->ops->view = PCView_FieldSplit_GKB;
2806: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);
2807: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetGKBTol_C",PCFieldSplitSetGKBTol_FieldSplit);
2808: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetGKBMaxit_C",PCFieldSplitSetGKBMaxit_FieldSplit);
2809: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetGKBNu_C",PCFieldSplitSetGKBNu_FieldSplit);
2810: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetGKBDelay_C",PCFieldSplitSetGKBDelay_FieldSplit);
2811: } else {
2812: pc->ops->apply = PCApply_FieldSplit;
2813: pc->ops->view = PCView_FieldSplit;
2815: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);
2816: }
2817: return(0);
2818: }
2820: static PetscErrorCode PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
2821: {
2822: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2825: if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
2826: if (jac->bs > 0 && jac->bs != bs) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change fieldsplit blocksize from %D to %D after it has been set",jac->bs,bs);
2827: jac->bs = bs;
2828: return(0);
2829: }
2831: /*@
2832: PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
2834: Collective on PC
2836: Input Parameter:
2837: + pc - the preconditioner context
2838: - type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
2840: Options Database Key:
2841: . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
2843: Level: Intermediate
2845: .seealso: PCCompositeSetType()
2847: @*/
2848: PetscErrorCode PCFieldSplitSetType(PC pc,PCCompositeType type)
2849: {
2854: PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));
2855: return(0);
2856: }
2858: /*@
2859: PCFieldSplitGetType - Gets the type of fieldsplit preconditioner.
2861: Not collective
2863: Input Parameter:
2864: . pc - the preconditioner context
2866: Output Parameter:
2867: . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
2869: Level: Intermediate
2871: .seealso: PCCompositeSetType()
2872: @*/
2873: PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type)
2874: {
2875: PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;
2880: *type = jac->type;
2881: return(0);
2882: }
2884: /*@
2885: PCFieldSplitSetDMSplits - Flags whether DMCreateFieldDecomposition() should be used to define the splits, whenever possible.
2887: Logically Collective
2889: Input Parameters:
2890: + pc - the preconditioner context
2891: - flg - boolean indicating whether to use field splits defined by the DM
2893: Options Database Key:
2894: . -pc_fieldsplit_dm_splits
2896: Level: Intermediate
2898: .seealso: PCFieldSplitGetDMSplits()
2900: @*/
2901: PetscErrorCode PCFieldSplitSetDMSplits(PC pc,PetscBool flg)
2902: {
2903: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2904: PetscBool isfs;
2910: PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
2911: if (isfs) {
2912: jac->dm_splits = flg;
2913: }
2914: return(0);
2915: }
2918: /*@
2919: PCFieldSplitGetDMSplits - Returns flag indicating whether DMCreateFieldDecomposition() should be used to define the splits, whenever possible.
2921: Logically Collective
2923: Input Parameter:
2924: . pc - the preconditioner context
2926: Output Parameter:
2927: . flg - boolean indicating whether to use field splits defined by the DM
2929: Level: Intermediate
2931: .seealso: PCFieldSplitSetDMSplits()
2933: @*/
2934: PetscErrorCode PCFieldSplitGetDMSplits(PC pc,PetscBool* flg)
2935: {
2936: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2937: PetscBool isfs;
2943: PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);
2944: if (isfs) {
2945: if(flg) *flg = jac->dm_splits;
2946: }
2947: return(0);
2948: }
2950: /*@
2951: PCFieldSplitGetDetectSaddlePoint - Returns flag indicating whether PCFieldSplit will attempt to automatically determine fields based on zero diagonal entries.
2953: Logically Collective
2955: Input Parameter:
2956: . pc - the preconditioner context
2958: Output Parameter:
2959: . flg - boolean indicating whether to detect fields or not
2961: Level: Intermediate
2963: .seealso: PCFIELDSPLIT, PCFieldSplitSetDetectSaddlePoint()
2965: @*/
2966: PetscErrorCode PCFieldSplitGetDetectSaddlePoint(PC pc,PetscBool *flg)
2967: {
2968: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2971: *flg = jac->detect;
2972: return(0);
2973: }
2975: /*@
2976: PCFieldSplitSetDetectSaddlePoint - Sets flag indicating whether PCFieldSplit will attempt to automatically determine fields based on zero diagonal entries.
2978: Logically Collective
2980: Notes:
2981: Also sets the split type to PC_COMPOSITE_SCHUR (see PCFieldSplitSetType()) and the Schur preconditioner type to PC_FIELDSPLIT_SCHUR_PRE_SELF (see PCFieldSplitSetSchurPre()).
2983: Input Parameter:
2984: . pc - the preconditioner context
2986: Output Parameter:
2987: . flg - boolean indicating whether to detect fields or not
2989: Options Database Key:
2990: . -pc_fieldsplit_detect_saddle_point
2992: Level: Intermediate
2994: .seealso: PCFIELDSPLIT, PCFieldSplitSetDetectSaddlePoint(), PCFieldSplitSetType(), PCFieldSplitSetSchurPre()
2996: @*/
2997: PetscErrorCode PCFieldSplitSetDetectSaddlePoint(PC pc,PetscBool flg)
2998: {
2999: PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
3003: jac->detect = flg;
3004: if (jac->detect) {
3005: PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);
3006: PCFieldSplitSetSchurPre(pc,PC_FIELDSPLIT_SCHUR_PRE_SELF,NULL);
3007: }
3008: return(0);
3009: }
3011: /* -------------------------------------------------------------------------------------*/
3012: /*MC
3013: PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
3014: fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.
3016: To set options on the solvers for each block append -fieldsplit_ to all the PC
3017: options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
3019: To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
3020: and set the options directly on the resulting KSP object
3022: Level: intermediate
3024: Options Database Keys:
3025: + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
3026: . -pc_fieldsplit_default - automatically add any fields to additional splits that have not
3027: been supplied explicitly by -pc_fieldsplit_%d_fields
3028: . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
3029: . -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur,gkb> - type of relaxation or factorization splitting
3030: . -pc_fieldsplit_schur_precondition <self,selfp,user,a11,full> - default is a11; see PCFieldSplitSetSchurPre()
3031: . -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> - set factorization type when using -pc_fieldsplit_type schur
3032: - -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero diagonal and uses Schur complement with no preconditioner as the solver
3034: Options prefixes for inner solvers when using the Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ .
3035: For all other solvers they are -fieldsplit_%d_ for the dth field; use -fieldsplit_ for all fields.
3036: The options prefix for the inner solver when using the Golub-Kahan biadiagonalization preconditioner is -fieldsplit_0_
3038: Notes:
3039: Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
3040: to define a field by an arbitrary collection of entries.
3042: If no fields are set the default is used. The fields are defined by entries strided by bs,
3043: beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
3044: if this is not called the block size defaults to the blocksize of the second matrix passed
3045: to KSPSetOperators()/PCSetOperators().
3047: $ For the Schur complement preconditioner if J = ( A00 A01 )
3048: $ ( A10 A11 )
3049: $ the preconditioner using full factorization is
3050: $ ( I -ksp(A00) A01 ) ( inv(A00) 0 ) ( I 0 )
3051: $ ( 0 I ) ( 0 ksp(S) ) ( -A10 ksp(A00) I )
3052: where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. S is the Schur complement
3053: $ S = A11 - A10 ksp(A00) A01
3054: which is usually dense and not stored explicitly. The action of ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given
3055: in providing the SECOND split or 1 if not give). For PCFieldSplitGetSubKSP() when field number is 0,
3056: it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
3057: A11 is used to construct a preconditioner for S, use PCFieldSplitSetSchurPre() for all the possible ways to construct the preconditioner for S.
3059: The factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above,
3060: diag gives
3061: $ ( inv(A00) 0 )
3062: $ ( 0 -ksp(S) )
3063: note that slightly counter intuitively there is a negative in front of the ksp(S) so that the preconditioner is positive definite. For SPD matrices J, the sign flip
3064: can be turned off with PCFieldSplitSetSchurScale() or by command line -pc_fieldsplit_schur_scale 1.0. The lower factorization is the inverse of
3065: $ ( A00 0 )
3066: $ ( A10 S )
3067: where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of
3068: $ ( A00 A01 )
3069: $ ( 0 S )
3070: where again the inverses of A00 and S are applied using KSPs.
3072: If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
3073: is used automatically for a second block.
3075: The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1.
3076: Generally it should be used with the AIJ format.
3078: The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see,
3079: for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT
3080: inside a smoother resulting in "Distributive Smoothers".
3082: There is a nice discussion of block preconditioners in
3084: [El08] A taxonomy and comparison of parallel block multi-level preconditioners for the incompressible Navier-Stokes equations
3085: Howard Elman, V.E. Howle, John Shadid, Robert Shuttleworth, Ray Tuminaro, Journal of Computational Physics 227 (2008) 1790--1808
3086: http://chess.cs.umd.edu/~elman/papers/tax.pdf
3088: The Constrained Pressure Preconditioner (CPR) can be implemented using PCCOMPOSITE with PCGALERKIN. CPR first solves an R A P subsystem, updates the
3089: residual on all variables (PCCompositeSetType(pc,PC_COMPOSITE_MULTIPLICATIVE)), and then applies a simple ILU like preconditioner on all the variables.
3091: The generalized Golub-Kahan bidiagonalization preconditioner (gkb) can be applied to symmetric 2x2 block matrices of the shape
3092: $ ( A00 A01 )
3093: $ ( A01' 0 )
3094: with A00 positive semi-definite. The implementation follows [Ar13]. Therein, we choose N := 1/nu * I and the (1,1)-block of the matrix is modified to H = A00 + nu*A01*A01'.
3095: A linear system Hx = b has to be solved in each iteration of the GKB algorithm. This solver is chosen with the option prefix -fieldsplit_0_.
3097: [Ar13] Generalized Golub-Kahan bidiagonalization and stopping criteria, SIAM J. Matrix Anal. Appl., Vol. 34, No. 2, pp. 571-592, 2013.
3099: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
3100: PCFieldSplitGetSubKSP(), PCFieldSplitSchurGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSetSchurPre(),
3101: MatSchurComplementSetAinvType(), PCFieldSplitSetSchurScale(),
3102: PCFieldSplitSetDetectSaddlePoint()
3103: M*/
3105: PETSC_EXTERN PetscErrorCode PCCreate_FieldSplit(PC pc)
3106: {
3108: PC_FieldSplit *jac;
3111: PetscNewLog(pc,&jac);
3113: jac->bs = -1;
3114: jac->nsplits = 0;
3115: jac->type = PC_COMPOSITE_MULTIPLICATIVE;
3116: jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
3117: jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL;
3118: jac->schurscale = -1.0;
3119: jac->dm_splits = PETSC_TRUE;
3120: jac->detect = PETSC_FALSE;
3121: jac->gkbtol = 1e-5;
3122: jac->gkbdelay = 5;
3123: jac->gkbnu = 1;
3124: jac->gkbmaxit = 100;
3125: jac->gkbmonitor = PETSC_FALSE;
3127: pc->data = (void*)jac;
3129: pc->ops->apply = PCApply_FieldSplit;
3130: pc->ops->applytranspose = PCApplyTranspose_FieldSplit;
3131: pc->ops->setup = PCSetUp_FieldSplit;
3132: pc->ops->reset = PCReset_FieldSplit;
3133: pc->ops->destroy = PCDestroy_FieldSplit;
3134: pc->ops->setfromoptions = PCSetFromOptions_FieldSplit;
3135: pc->ops->view = PCView_FieldSplit;
3136: pc->ops->applyrichardson = 0;
3138: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSchurGetSubKSP_C",PCFieldSplitSchurGetSubKSP_FieldSplit);
3139: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);
3140: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",PCFieldSplitSetFields_FieldSplit);
3141: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",PCFieldSplitSetIS_FieldSplit);
3142: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",PCFieldSplitSetType_FieldSplit);
3143: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",PCFieldSplitSetBlockSize_FieldSplit);
3144: PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitRestrictIS_C",PCFieldSplitRestrictIS_FieldSplit);
3145: return(0);
3146: }