xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision a80b646ebe2cda77ae9dde47e1d7b0a0f620364d)
1dba47a55SKris Buschelman 
2af0996ceSBarry Smith #include <petsc/private/pcimpl.h>     /*I "petscpc.h" I*/
3*a80b646eSBarry Smith #include <petsc/private/kspimpl.h>    /*  This is needed to provide the appropriate PETSC_EXTERN for KSP_Solve_FS ....*/
41e25c274SJed Brown #include <petscdm.h>
50971522cSBarry Smith 
62e71c61dSMatthew G. Knepley const char *const PCFieldSplitSchurPreTypes[] = {"SELF","SELFP","A11","USER","FULL","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0};
7c9c6ffaaSJed Brown const char *const PCFieldSplitSchurFactTypes[] = {"DIAG","LOWER","UPPER","FULL","PCFieldSplitSchurFactType","PC_FIELDSPLIT_SCHUR_FACT_",0};
8c5d2311dSJed Brown 
99df09d43SBarry Smith 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;
109df09d43SBarry Smith 
110971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
120971522cSBarry Smith struct _PC_FieldSplitLink {
1369a612a9SBarry Smith   KSP               ksp;
14443836d0SMatthew G Knepley   Vec               x,y,z;
15db4c96c1SJed Brown   char              *splitname;
160971522cSBarry Smith   PetscInt          nfields;
175d4c12cdSJungho Lee   PetscInt          *fields,*fields_col;
181b9fc7fcSBarry Smith   VecScatter        sctx;
196dbb499eSCian Wilson   IS                is,is_col,is_orig;
2051f519a2SBarry Smith   PC_FieldSplitLink next,previous;
219df09d43SBarry Smith   PetscLogEvent     event;
220971522cSBarry Smith };
230971522cSBarry Smith 
240971522cSBarry Smith typedef struct {
2568dd23aaSBarry Smith   PCCompositeType type;
26ace3abfcSBarry Smith   PetscBool       defaultsplit;                    /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */
27ace3abfcSBarry Smith   PetscBool       splitdefined;                    /* Flag is set after the splits have been defined, to prevent more splits from being added */
2830ad9308SMatthew Knepley   PetscInt        bs;                              /* Block size for IS and Mat structures */
2930ad9308SMatthew Knepley   PetscInt        nsplits;                         /* Number of field divisions defined */
3079416396SBarry Smith   Vec             *x,*y,w1,w2;
31519d70e2SJed Brown   Mat             *mat;                            /* The diagonal block for each split */
32519d70e2SJed Brown   Mat             *pmat;                           /* The preconditioning diagonal block for each split */
3330ad9308SMatthew Knepley   Mat             *Afield;                         /* The rows of the matrix associated with each split */
34ace3abfcSBarry Smith   PetscBool       issetup;
352fa5cd67SKarl Rupp 
3630ad9308SMatthew Knepley   /* Only used when Schur complement preconditioning is used */
3730ad9308SMatthew Knepley   Mat                       B;                     /* The (0,1) block */
3830ad9308SMatthew Knepley   Mat                       C;                     /* The (1,0) block */
39443836d0SMatthew G Knepley   Mat                       schur;                 /* The Schur complement S = A11 - A10 A00^{-1} A01, the KSP here, kspinner, is H_1 in [El08] */
40a7476a74SDmitry Karpeev   Mat                       schurp;                /* Assembled approximation to S built by MatSchurComplement to be used as a preconditioning matrix when solving with S */
41084e4875SJed Brown   Mat                       schur_user;            /* User-provided preconditioning matrix for the Schur complement */
42084e4875SJed Brown   PCFieldSplitSchurPreType  schurpre;              /* Determines which preconditioning matrix is used for the Schur complement */
43c9c6ffaaSJed Brown   PCFieldSplitSchurFactType schurfactorization;
4430ad9308SMatthew Knepley   KSP                       kspschur;              /* The solver for S */
45443836d0SMatthew G Knepley   KSP                       kspupper;              /* The solver for A in the upper diagonal part of the factorization (H_2 in [El08]) */
4697bbdb24SBarry Smith   PC_FieldSplitLink         head;
4763ec74ffSBarry Smith   PetscBool                 reset;                  /* indicates PCReset() has been last called on this object, hack */
486dbb499eSCian Wilson   PetscBool                 isrestrict;             /* indicates PCFieldSplitRestrictIS() has been last called on this object, hack */
49c1570756SJed Brown   PetscBool                 suboptionsset;          /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */
504ab8060aSDmitry Karpeev   PetscBool                 dm_splits;              /* Whether to use DMCreateFieldDecomposition() whenever possible */
51c84da90fSDmitry Karpeev   PetscBool                 diag_use_amat;          /* Whether to extract diagonal matrix blocks from Amat, rather than Pmat (weaker than -pc_use_amat) */
52c84da90fSDmitry Karpeev   PetscBool                 offdiag_use_amat;       /* Whether to extract off-diagonal matrix blocks from Amat, rather than Pmat (weaker than -pc_use_amat) */
530971522cSBarry Smith } PC_FieldSplit;
540971522cSBarry Smith 
5516913363SBarry Smith /*
5616913363SBarry Smith     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
5716913363SBarry Smith    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
5816913363SBarry Smith    PC you could change this.
5916913363SBarry Smith */
60084e4875SJed Brown 
61e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it.  This way the
62084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */
63084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac)
64084e4875SJed Brown {
65084e4875SJed Brown   switch (jac->schurpre) {
66084e4875SJed Brown   case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur;
67a7476a74SDmitry Karpeev   case PC_FIELDSPLIT_SCHUR_PRE_SELFP: return jac->schurp;
68e87fab1bSBarry Smith   case PC_FIELDSPLIT_SCHUR_PRE_A11: return jac->pmat[1];
69e74569cdSMatthew G. Knepley   case PC_FIELDSPLIT_SCHUR_PRE_FULL: /* We calculate this and store it in schur_user */
70084e4875SJed Brown   case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */
71084e4875SJed Brown   default:
72084e4875SJed Brown     return jac->schur_user ? jac->schur_user : jac->pmat[1];
73084e4875SJed Brown   }
74084e4875SJed Brown }
75084e4875SJed Brown 
76084e4875SJed Brown 
779804daf3SBarry Smith #include <petscdraw.h>
780971522cSBarry Smith #undef __FUNCT__
790971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
800971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
810971522cSBarry Smith {
820971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
830971522cSBarry Smith   PetscErrorCode    ierr;
84d9884438SBarry Smith   PetscBool         iascii,isdraw;
850971522cSBarry Smith   PetscInt          i,j;
865a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
870971522cSBarry Smith 
880971522cSBarry Smith   PetscFunctionBegin;
89251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
90d9884438SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
910971522cSBarry Smith   if (iascii) {
922c9966d7SBarry Smith     if (jac->bs > 0) {
9351f519a2SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
942c9966d7SBarry Smith     } else {
952c9966d7SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D\n",PCCompositeTypes[jac->type],jac->nsplits);CHKERRQ(ierr);
962c9966d7SBarry Smith     }
97f5236f50SJed Brown     if (pc->useAmat) {
98f5236f50SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for blocks\n");CHKERRQ(ierr);
99a3df900dSMatthew G Knepley     }
100c84da90fSDmitry Karpeev     if (jac->diag_use_amat) {
101c84da90fSDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for diagonal blocks\n");CHKERRQ(ierr);
102c84da90fSDmitry Karpeev     }
103c84da90fSDmitry Karpeev     if (jac->offdiag_use_amat) {
104c84da90fSDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for off-diagonal blocks\n");CHKERRQ(ierr);
105c84da90fSDmitry Karpeev     }
10669a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
1070971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1080971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
1091ab39975SBarry Smith       if (ilink->fields) {
1100971522cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
11179416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
1125a9f2f41SSatish Balay         for (j=0; j<ilink->nfields; j++) {
11379416396SBarry Smith           if (j > 0) {
11479416396SBarry Smith             ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
11579416396SBarry Smith           }
1165a9f2f41SSatish Balay           ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
1170971522cSBarry Smith         }
1180971522cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
11979416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1201ab39975SBarry Smith       } else {
1211ab39975SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1221ab39975SBarry Smith       }
1235a9f2f41SSatish Balay       ierr  = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
1245a9f2f41SSatish Balay       ilink = ilink->next;
1250971522cSBarry Smith     }
1260971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1272fa5cd67SKarl Rupp   }
1282fa5cd67SKarl Rupp 
1292fa5cd67SKarl Rupp  if (isdraw) {
130d9884438SBarry Smith     PetscDraw draw;
131d9884438SBarry Smith     PetscReal x,y,w,wd;
132d9884438SBarry Smith 
133d9884438SBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
134d9884438SBarry Smith     ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
135d9884438SBarry Smith     w    = 2*PetscMin(1.0 - x,x);
136d9884438SBarry Smith     wd   = w/(jac->nsplits + 1);
137d9884438SBarry Smith     x    = x - wd*(jac->nsplits-1)/2.0;
138d9884438SBarry Smith     for (i=0; i<jac->nsplits; i++) {
139d9884438SBarry Smith       ierr  = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
140d9884438SBarry Smith       ierr  = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
141d9884438SBarry Smith       ierr  = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
142d9884438SBarry Smith       x    += wd;
143d9884438SBarry Smith       ilink = ilink->next;
144d9884438SBarry Smith     }
1450971522cSBarry Smith   }
1460971522cSBarry Smith   PetscFunctionReturn(0);
1470971522cSBarry Smith }
1480971522cSBarry Smith 
1490971522cSBarry Smith #undef __FUNCT__
1503b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur"
1513b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
1523b224e63SBarry Smith {
1533b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1543b224e63SBarry Smith   PetscErrorCode    ierr;
1554996c5bdSBarry Smith   PetscBool         iascii,isdraw;
1563b224e63SBarry Smith   PetscInt          i,j;
1573b224e63SBarry Smith   PC_FieldSplitLink ilink = jac->head;
1583b224e63SBarry Smith 
1593b224e63SBarry Smith   PetscFunctionBegin;
160251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1614996c5bdSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
1623b224e63SBarry Smith   if (iascii) {
1632c9966d7SBarry Smith     if (jac->bs > 0) {
164c9c6ffaaSJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
1652c9966d7SBarry Smith     } else {
1662c9966d7SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, factorization %s\n",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
1672c9966d7SBarry Smith     }
168f5236f50SJed Brown     if (pc->useAmat) {
169f5236f50SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for blocks\n");CHKERRQ(ierr);
170a3df900dSMatthew G Knepley     }
1713e8b8b31SMatthew G Knepley     switch (jac->schurpre) {
1723e8b8b31SMatthew G Knepley     case PC_FIELDSPLIT_SCHUR_PRE_SELF:
1733e8b8b31SMatthew G Knepley       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from S itself\n");CHKERRQ(ierr);break;
174a7476a74SDmitry Karpeev     case PC_FIELDSPLIT_SCHUR_PRE_SELFP:
1751d71051fSDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from Sp, an assembled approximation to S, which uses (lumped, if requested) A00's diagonal's inverse\n");CHKERRQ(ierr);break;
176e87fab1bSBarry Smith     case PC_FIELDSPLIT_SCHUR_PRE_A11:
177e87fab1bSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from A11\n");CHKERRQ(ierr);break;
178e74569cdSMatthew G. Knepley     case PC_FIELDSPLIT_SCHUR_PRE_FULL:
179e74569cdSMatthew G. Knepley       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from the exact Schur complement\n");CHKERRQ(ierr);break;
1803e8b8b31SMatthew G Knepley     case PC_FIELDSPLIT_SCHUR_PRE_USER:
1813e8b8b31SMatthew G Knepley       if (jac->schur_user) {
1823e8b8b31SMatthew G Knepley         ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from user provided matrix\n");CHKERRQ(ierr);
1833e8b8b31SMatthew G Knepley       } else {
184e87fab1bSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from A11\n");CHKERRQ(ierr);
1853e8b8b31SMatthew G Knepley       }
1863e8b8b31SMatthew G Knepley       break;
1873e8b8b31SMatthew G Knepley     default:
18882f516ccSBarry Smith       SETERRQ1(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Invalid Schur preconditioning type: %d", jac->schurpre);
1893e8b8b31SMatthew G Knepley     }
1903b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Split info:\n");CHKERRQ(ierr);
1913b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1923b224e63SBarry Smith     for (i=0; i<jac->nsplits; i++) {
1933b224e63SBarry Smith       if (ilink->fields) {
1943b224e63SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
1953b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
1963b224e63SBarry Smith         for (j=0; j<ilink->nfields; j++) {
1973b224e63SBarry Smith           if (j > 0) {
1983b224e63SBarry Smith             ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
1993b224e63SBarry Smith           }
2003b224e63SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
2013b224e63SBarry Smith         }
2023b224e63SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
2033b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
2043b224e63SBarry Smith       } else {
2053b224e63SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
2063b224e63SBarry Smith       }
2073b224e63SBarry Smith       ilink = ilink->next;
2083b224e63SBarry Smith     }
209435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block\n");CHKERRQ(ierr);
2103b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
21106de4afeSJed Brown     if (jac->head) {
212443836d0SMatthew G Knepley       ierr = KSPView(jac->head->ksp,viewer);CHKERRQ(ierr);
21306de4afeSJed Brown     } else  {ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);}
2143b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
21506de4afeSJed Brown     if (jac->head && jac->kspupper != jac->head->ksp) {
216443836d0SMatthew G Knepley       ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for upper A00 in upper triangular factor \n");CHKERRQ(ierr);
217443836d0SMatthew G Knepley       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
218b2750c55SJed Brown       if (jac->kspupper) {ierr = KSPView(jac->kspupper,viewer);CHKERRQ(ierr);}
219b2750c55SJed Brown       else {ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);}
220443836d0SMatthew G Knepley       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
221443836d0SMatthew G Knepley     }
222435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr);
2233b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
22412cae6f2SJed Brown     if (jac->kspschur) {
2253b224e63SBarry Smith       ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
22612cae6f2SJed Brown     } else {
22712cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
22812cae6f2SJed Brown     }
2293b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
2303b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
23106de4afeSJed Brown   } else if (isdraw && jac->head) {
2324996c5bdSBarry Smith     PetscDraw draw;
2334996c5bdSBarry Smith     PetscReal x,y,w,wd,h;
2344996c5bdSBarry Smith     PetscInt  cnt = 2;
2354996c5bdSBarry Smith     char      str[32];
2364996c5bdSBarry Smith 
2374996c5bdSBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
2384996c5bdSBarry Smith     ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
239c74581afSBarry Smith     if (jac->kspupper != jac->head->ksp) cnt++;
240c74581afSBarry Smith     w  = 2*PetscMin(1.0 - x,x);
241c74581afSBarry Smith     wd = w/(cnt + 1);
242c74581afSBarry Smith 
2434996c5bdSBarry Smith     ierr = PetscSNPrintf(str,32,"Schur fact. %s",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
24451fa3d41SBarry Smith     ierr = PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr);
2454996c5bdSBarry Smith     y   -= h;
2464996c5bdSBarry Smith     if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_USER &&  !jac->schur_user) {
247e87fab1bSBarry Smith       ierr = PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[PC_FIELDSPLIT_SCHUR_PRE_A11]);CHKERRQ(ierr);
2483b224e63SBarry Smith     } else {
2494996c5bdSBarry Smith       ierr = PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[jac->schurpre]);CHKERRQ(ierr);
2504996c5bdSBarry Smith     }
25151fa3d41SBarry Smith     ierr = PetscDrawStringBoxed(draw,x+wd*(cnt-1)/2.0,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr);
2524996c5bdSBarry Smith     y   -= h;
2534996c5bdSBarry Smith     x    = x - wd*(cnt-1)/2.0;
2544996c5bdSBarry Smith 
2554996c5bdSBarry Smith     ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
2564996c5bdSBarry Smith     ierr = KSPView(jac->head->ksp,viewer);CHKERRQ(ierr);
2574996c5bdSBarry Smith     ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
2584996c5bdSBarry Smith     if (jac->kspupper != jac->head->ksp) {
2594996c5bdSBarry Smith       x   += wd;
2604996c5bdSBarry Smith       ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
2614996c5bdSBarry Smith       ierr = KSPView(jac->kspupper,viewer);CHKERRQ(ierr);
2624996c5bdSBarry Smith       ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
2634996c5bdSBarry Smith     }
2644996c5bdSBarry Smith     x   += wd;
2654996c5bdSBarry Smith     ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
2664996c5bdSBarry Smith     ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
2674996c5bdSBarry Smith     ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
2683b224e63SBarry Smith   }
2693b224e63SBarry Smith   PetscFunctionReturn(0);
2703b224e63SBarry Smith }
2713b224e63SBarry Smith 
2723b224e63SBarry Smith #undef __FUNCT__
2736c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private"
2746c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */
2756c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc)
2766c924f48SJed Brown {
2776c924f48SJed Brown   PetscErrorCode ierr;
2786c924f48SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
2795d4c12cdSJungho Lee   PetscInt       i,nfields,*ifields,nfields_col,*ifields_col;
2805d4c12cdSJungho Lee   PetscBool      flg,flg_col;
2815d4c12cdSJungho Lee   char           optionname[128],splitname[8],optionname_col[128];
2826c924f48SJed Brown 
2836c924f48SJed Brown   PetscFunctionBegin;
284785e854fSJed Brown   ierr = PetscMalloc1(jac->bs,&ifields);CHKERRQ(ierr);
285785e854fSJed Brown   ierr = PetscMalloc1(jac->bs,&ifields_col);CHKERRQ(ierr);
2866c924f48SJed Brown   for (i=0,flg=PETSC_TRUE;; i++) {
2878caf3d72SBarry Smith     ierr        = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr);
2888caf3d72SBarry Smith     ierr        = PetscSNPrintf(optionname,sizeof(optionname),"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr);
2898caf3d72SBarry Smith     ierr        = PetscSNPrintf(optionname_col,sizeof(optionname_col),"-pc_fieldsplit_%D_fields_col",i);CHKERRQ(ierr);
2906c924f48SJed Brown     nfields     = jac->bs;
29129499fbbSJungho Lee     nfields_col = jac->bs;
292c5929fdfSBarry Smith     ierr        = PetscOptionsGetIntArray(((PetscObject)pc)->options,((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr);
293c5929fdfSBarry Smith     ierr        = PetscOptionsGetIntArray(((PetscObject)pc)->options,((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);CHKERRQ(ierr);
2946c924f48SJed Brown     if (!flg) break;
2955d4c12cdSJungho Lee     else if (flg && !flg_col) {
2965d4c12cdSJungho Lee       if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
2975d4c12cdSJungho Lee       ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);CHKERRQ(ierr);
2982fa5cd67SKarl Rupp     } else {
2995d4c12cdSJungho Lee       if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
3005d4c12cdSJungho Lee       if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match");
3015d4c12cdSJungho Lee       ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);CHKERRQ(ierr);
3025d4c12cdSJungho Lee     }
3036c924f48SJed Brown   }
3046c924f48SJed Brown   if (i > 0) {
3056c924f48SJed Brown     /* Makes command-line setting of splits take precedence over setting them in code.
3066c924f48SJed Brown        Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would
3076c924f48SJed Brown        create new splits, which would probably not be what the user wanted. */
3086c924f48SJed Brown     jac->splitdefined = PETSC_TRUE;
3096c924f48SJed Brown   }
3106c924f48SJed Brown   ierr = PetscFree(ifields);CHKERRQ(ierr);
3115d4c12cdSJungho Lee   ierr = PetscFree(ifields_col);CHKERRQ(ierr);
3126c924f48SJed Brown   PetscFunctionReturn(0);
3136c924f48SJed Brown }
3146c924f48SJed Brown 
3156c924f48SJed Brown #undef __FUNCT__
31669a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
31769a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
3180971522cSBarry Smith {
3190971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
3200971522cSBarry Smith   PetscErrorCode    ierr;
3215a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
3223a062f41SBarry Smith   PetscBool         fieldsplit_default = PETSC_FALSE,stokes = PETSC_FALSE,coupling = PETSC_FALSE;
3236c924f48SJed Brown   PetscInt          i;
3240971522cSBarry Smith 
3250971522cSBarry Smith   PetscFunctionBegin;
3267287d2a3SDmitry Karpeev   /*
3277287d2a3SDmitry Karpeev    Kinda messy, but at least this now uses DMCreateFieldDecomposition() even with jac->reset.
3287287d2a3SDmitry Karpeev    Should probably be rewritten.
3297287d2a3SDmitry Karpeev    */
3307287d2a3SDmitry Karpeev   if (!ilink || jac->reset) {
331c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);CHKERRQ(ierr);
332c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_coupling",&coupling,NULL);CHKERRQ(ierr);
3333a062f41SBarry Smith     if (pc->dm && jac->dm_splits && !stokes && !coupling) {
334bafc1b83SMatthew G Knepley       PetscInt  numFields, f, i, j;
3350784a22cSJed Brown       char      **fieldNames;
3367b62db95SJungho Lee       IS        *fields;
337e7c4fc90SDmitry Karpeev       DM        *dms;
338bafc1b83SMatthew G Knepley       DM        subdm[128];
339bafc1b83SMatthew G Knepley       PetscBool flg;
340bafc1b83SMatthew G Knepley 
34116621825SDmitry Karpeev       ierr = DMCreateFieldDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms);CHKERRQ(ierr);
342bafc1b83SMatthew G Knepley       /* Allow the user to prescribe the splits */
343bafc1b83SMatthew G Knepley       for (i = 0, flg = PETSC_TRUE;; i++) {
344bafc1b83SMatthew G Knepley         PetscInt ifields[128];
345bafc1b83SMatthew G Knepley         IS       compField;
346bafc1b83SMatthew G Knepley         char     optionname[128], splitname[8];
347bafc1b83SMatthew G Knepley         PetscInt nfields = numFields;
348bafc1b83SMatthew G Knepley 
3498caf3d72SBarry Smith         ierr = PetscSNPrintf(optionname, sizeof(optionname), "-pc_fieldsplit_%D_fields", i);CHKERRQ(ierr);
350c5929fdfSBarry Smith         ierr = PetscOptionsGetIntArray(((PetscObject)pc)->options,((PetscObject)pc)->prefix, optionname, ifields, &nfields, &flg);CHKERRQ(ierr);
351bafc1b83SMatthew G Knepley         if (!flg) break;
35282f516ccSBarry Smith         if (numFields > 128) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot currently support %d > 128 fields", numFields);
353bafc1b83SMatthew G Knepley         ierr = DMCreateSubDM(pc->dm, nfields, ifields, &compField, &subdm[i]);CHKERRQ(ierr);
354bafc1b83SMatthew G Knepley         if (nfields == 1) {
355bafc1b83SMatthew G Knepley           ierr = PCFieldSplitSetIS(pc, fieldNames[ifields[0]], compField);CHKERRQ(ierr);
356bafc1b83SMatthew G Knepley         } else {
3578caf3d72SBarry Smith           ierr = PetscSNPrintf(splitname, sizeof(splitname), "%D", i);CHKERRQ(ierr);
358bafc1b83SMatthew G Knepley           ierr = PCFieldSplitSetIS(pc, splitname, compField);CHKERRQ(ierr);
3597287d2a3SDmitry Karpeev         }
360bafc1b83SMatthew G Knepley         ierr = ISDestroy(&compField);CHKERRQ(ierr);
361bafc1b83SMatthew G Knepley         for (j = 0; j < nfields; ++j) {
362bafc1b83SMatthew G Knepley           f    = ifields[j];
3637b62db95SJungho Lee           ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr);
3647b62db95SJungho Lee           ierr = ISDestroy(&fields[f]);CHKERRQ(ierr);
3657b62db95SJungho Lee         }
366bafc1b83SMatthew G Knepley       }
367bafc1b83SMatthew G Knepley       if (i == 0) {
368bafc1b83SMatthew G Knepley         for (f = 0; f < numFields; ++f) {
369bafc1b83SMatthew G Knepley           ierr = PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);CHKERRQ(ierr);
370bafc1b83SMatthew G Knepley           ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr);
371bafc1b83SMatthew G Knepley           ierr = ISDestroy(&fields[f]);CHKERRQ(ierr);
372bafc1b83SMatthew G Knepley         }
373bafc1b83SMatthew G Knepley       } else {
374d724dfffSBarry Smith         for (j=0; j<numFields; j++) {
375d724dfffSBarry Smith           ierr = DMDestroy(dms+j);CHKERRQ(ierr);
376d724dfffSBarry Smith         }
377d724dfffSBarry Smith         ierr = PetscFree(dms);CHKERRQ(ierr);
378785e854fSJed Brown         ierr = PetscMalloc1(i, &dms);CHKERRQ(ierr);
3792fa5cd67SKarl Rupp         for (j = 0; j < i; ++j) dms[j] = subdm[j];
380bafc1b83SMatthew G Knepley       }
3817b62db95SJungho Lee       ierr = PetscFree(fieldNames);CHKERRQ(ierr);
3827b62db95SJungho Lee       ierr = PetscFree(fields);CHKERRQ(ierr);
383e7c4fc90SDmitry Karpeev       if (dms) {
3848b8307b2SJed Brown         ierr = PetscInfo(pc, "Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
385bafc1b83SMatthew G Knepley         for (ilink = jac->head, i = 0; ilink; ilink = ilink->next, ++i) {
3867287d2a3SDmitry Karpeev           const char *prefix;
3877287d2a3SDmitry Karpeev           ierr = PetscObjectGetOptionsPrefix((PetscObject)(ilink->ksp),&prefix);CHKERRQ(ierr);
3887287d2a3SDmitry Karpeev           ierr = PetscObjectSetOptionsPrefix((PetscObject)(dms[i]), prefix);CHKERRQ(ierr);
3897b62db95SJungho Lee           ierr = KSPSetDM(ilink->ksp, dms[i]);CHKERRQ(ierr);
3907b62db95SJungho Lee           ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE);CHKERRQ(ierr);
3917287d2a3SDmitry Karpeev           ierr = PetscObjectIncrementTabLevel((PetscObject)dms[i],(PetscObject)ilink->ksp,0);CHKERRQ(ierr);
392e7c4fc90SDmitry Karpeev           ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);
3932fa5ba8aSJed Brown         }
3947b62db95SJungho Lee         ierr = PetscFree(dms);CHKERRQ(ierr);
3958b8307b2SJed Brown       }
39666ffff09SJed Brown     } else {
397521d7252SBarry Smith       if (jac->bs <= 0) {
398704ba839SBarry Smith         if (pc->pmat) {
399521d7252SBarry Smith           ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
4002fa5cd67SKarl Rupp         } else jac->bs = 1;
401521d7252SBarry Smith       }
402d32f9abdSBarry Smith 
4036ce1633cSBarry Smith       if (stokes) {
4046ce1633cSBarry Smith         IS       zerodiags,rest;
4056ce1633cSBarry Smith         PetscInt nmin,nmax;
4066ce1633cSBarry Smith 
4076ce1633cSBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
4086ce1633cSBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
4096ce1633cSBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
4107287d2a3SDmitry Karpeev         if (jac->reset) {
4112a808120SBarry Smith           if (!jac->head) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"At reset jac must have head from previous setup");
4127287d2a3SDmitry Karpeev           jac->head->is       = rest;
4137287d2a3SDmitry Karpeev           jac->head->next->is = zerodiags;
4142fa5cd67SKarl Rupp         } else {
4156ce1633cSBarry Smith           ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
4166ce1633cSBarry Smith           ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr);
4177287d2a3SDmitry Karpeev         }
418fcfd50ebSBarry Smith         ierr = ISDestroy(&zerodiags);CHKERRQ(ierr);
419fcfd50ebSBarry Smith         ierr = ISDestroy(&rest);CHKERRQ(ierr);
4203a062f41SBarry Smith       } else if (coupling) {
4213a062f41SBarry Smith         IS       coupling,rest;
4223a062f41SBarry Smith         PetscInt nmin,nmax;
4233a062f41SBarry Smith 
4243a062f41SBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
4253a062f41SBarry Smith         ierr = MatFindOffBlockDiagonalEntries(pc->mat,&coupling);CHKERRQ(ierr);
4266bea0878SBarry Smith         ierr = ISCreateStride(PetscObjectComm((PetscObject)pc->mat),nmax-nmin,nmin,1,&rest);CHKERRQ(ierr);
427e52d2c62SBarry Smith         ierr = ISSetIdentity(rest);CHKERRQ(ierr);
4283a062f41SBarry Smith         if (jac->reset) {
4292a808120SBarry Smith           if (!jac->head) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"At reset jac must have head from previous setup");
430bb3d4b4cSBarry Smith           jac->head->is       = rest;
431bb3d4b4cSBarry Smith           jac->head->next->is = coupling;
4323a062f41SBarry Smith         } else {
433d020c841SBarry Smith           ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
434d020c841SBarry Smith           ierr = PCFieldSplitSetIS(pc,"1",coupling);CHKERRQ(ierr);
4353a062f41SBarry Smith         }
4363a062f41SBarry Smith         ierr = ISDestroy(&coupling);CHKERRQ(ierr);
4373a062f41SBarry Smith         ierr = ISDestroy(&rest);CHKERRQ(ierr);
4386ce1633cSBarry Smith       } else {
4396dbb499eSCian Wilson         if (jac->reset && !jac->isrestrict) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used");
440c5929fdfSBarry Smith         ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&fieldsplit_default,NULL);CHKERRQ(ierr);
4417287d2a3SDmitry Karpeev         if (!fieldsplit_default) {
442d32f9abdSBarry Smith           /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
443d32f9abdSBarry Smith            then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
4446c924f48SJed Brown           ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
4456c924f48SJed Brown           if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
446d32f9abdSBarry Smith         }
4476dbb499eSCian Wilson         if ((fieldsplit_default || !jac->splitdefined) && !jac->isrestrict) {
448d32f9abdSBarry Smith           ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
449db4c96c1SJed Brown           for (i=0; i<jac->bs; i++) {
4506c924f48SJed Brown             char splitname[8];
4518caf3d72SBarry Smith             ierr = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr);
4525d4c12cdSJungho Lee             ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr);
45379416396SBarry Smith           }
4545d4c12cdSJungho Lee           jac->defaultsplit = PETSC_TRUE;
455521d7252SBarry Smith         }
45666ffff09SJed Brown       }
4576ce1633cSBarry Smith     }
458edf189efSBarry Smith   } else if (jac->nsplits == 1) {
459edf189efSBarry Smith     if (ilink->is) {
460edf189efSBarry Smith       IS       is2;
461edf189efSBarry Smith       PetscInt nmin,nmax;
462edf189efSBarry Smith 
463edf189efSBarry Smith       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
464edf189efSBarry Smith       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
465db4c96c1SJed Brown       ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr);
466fcfd50ebSBarry Smith       ierr = ISDestroy(&is2);CHKERRQ(ierr);
46782f516ccSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
468edf189efSBarry Smith   }
469d0af7cd3SBarry Smith 
47082f516ccSBarry Smith   if (jac->nsplits < 2) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits);
47169a612a9SBarry Smith   PetscFunctionReturn(0);
47269a612a9SBarry Smith }
47369a612a9SBarry Smith 
474c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(PetscOptions,const char pre[], const char name[], char *value[], PetscBool *flg);
475514bf10dSMatthew G Knepley 
47669a612a9SBarry Smith #undef __FUNCT__
47769a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
47869a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
47969a612a9SBarry Smith {
48069a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
48169a612a9SBarry Smith   PetscErrorCode    ierr;
4825a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
4832c9966d7SBarry Smith   PetscInt          i,nsplit;
4845f4ab4e1SJungho Lee   PetscBool         sorted, sorted_col;
48569a612a9SBarry Smith 
48669a612a9SBarry Smith   PetscFunctionBegin;
48769a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
48897bbdb24SBarry Smith   nsplit = jac->nsplits;
4895a9f2f41SSatish Balay   ilink  = jac->head;
49097bbdb24SBarry Smith 
49197bbdb24SBarry Smith   /* get the matrices for each split */
492704ba839SBarry Smith   if (!jac->issetup) {
4931b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
49497bbdb24SBarry Smith 
495704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
496704ba839SBarry Smith 
4975d4c12cdSJungho Lee     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
4982c9966d7SBarry Smith     if (jac->defaultsplit || !ilink->is) {
4992c9966d7SBarry Smith       if (jac->bs <= 0) jac->bs = nsplit;
5002c9966d7SBarry Smith     }
50151f519a2SBarry Smith     bs     = jac->bs;
50297bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
5031b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
5041b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
5051b9fc7fcSBarry Smith       if (jac->defaultsplit) {
506ce94432eSBarry Smith         ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
5075f4ab4e1SJungho Lee         ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr);
508704ba839SBarry Smith       } else if (!ilink->is) {
509ccb205f8SBarry Smith         if (ilink->nfields > 1) {
5105f4ab4e1SJungho Lee           PetscInt *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col;
511785e854fSJed Brown           ierr = PetscMalloc1(ilink->nfields*nslots,&ii);CHKERRQ(ierr);
512785e854fSJed Brown           ierr = PetscMalloc1(ilink->nfields*nslots,&jj);CHKERRQ(ierr);
5131b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
5141b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
5151b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
5165f4ab4e1SJungho Lee               jj[nfields*j + k] = rstart + bs*j + fields_col[k];
51797bbdb24SBarry Smith             }
51897bbdb24SBarry Smith           }
519ce94432eSBarry Smith           ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr);
520ce94432eSBarry Smith           ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr);
52190e68f20SPatrick Farrell           ierr = ISSetBlockSize(ilink->is, nfields);CHKERRQ(ierr);
52290e68f20SPatrick Farrell           ierr = ISSetBlockSize(ilink->is_col, nfields);CHKERRQ(ierr);
523ccb205f8SBarry Smith         } else {
524ce94432eSBarry Smith           ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
525ce94432eSBarry Smith           ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr);
526ccb205f8SBarry Smith        }
5273e197d65SBarry Smith       }
528704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
5295f4ab4e1SJungho Lee       if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); }
5305f4ab4e1SJungho Lee       if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
531704ba839SBarry Smith       ilink = ilink->next;
5321b9fc7fcSBarry Smith     }
5331b9fc7fcSBarry Smith   }
5341b9fc7fcSBarry Smith 
535704ba839SBarry Smith   ilink = jac->head;
53697bbdb24SBarry Smith   if (!jac->pmat) {
537bdddcaaaSMatthew G Knepley     Vec xtmp;
538bdddcaaaSMatthew G Knepley 
5392a7a6963SBarry Smith     ierr = MatCreateVecs(pc->pmat,&xtmp,NULL);CHKERRQ(ierr);
540785e854fSJed Brown     ierr = PetscMalloc1(nsplit,&jac->pmat);CHKERRQ(ierr);
541dcca6d9dSJed Brown     ierr = PetscMalloc2(nsplit,&jac->x,nsplit,&jac->y);CHKERRQ(ierr);
542cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
543bdddcaaaSMatthew G Knepley       MatNullSpace sp;
544bdddcaaaSMatthew G Knepley 
545a3df900dSMatthew G Knepley       /* Check for preconditioning matrix attached to IS */
546a3df900dSMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &jac->pmat[i]);CHKERRQ(ierr);
547a3df900dSMatthew G Knepley       if (jac->pmat[i]) {
548a3df900dSMatthew G Knepley         ierr = PetscObjectReference((PetscObject) jac->pmat[i]);CHKERRQ(ierr);
549a3df900dSMatthew G Knepley         if (jac->type == PC_COMPOSITE_SCHUR) {
550a3df900dSMatthew G Knepley           jac->schur_user = jac->pmat[i];
5512fa5cd67SKarl Rupp 
552a3df900dSMatthew G Knepley           ierr = PetscObjectReference((PetscObject) jac->schur_user);CHKERRQ(ierr);
553a3df900dSMatthew G Knepley         }
554a3df900dSMatthew G Knepley       } else {
5553a062f41SBarry Smith         const char *prefix;
5565f4ab4e1SJungho Lee         ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
5573a062f41SBarry Smith         ierr = KSPGetOptionsPrefix(ilink->ksp,&prefix);CHKERRQ(ierr);
5583a062f41SBarry Smith         ierr = MatSetOptionsPrefix(jac->pmat[i],prefix);CHKERRQ(ierr);
5593a062f41SBarry Smith         ierr = MatViewFromOptions(jac->pmat[i],NULL,"-mat_view");CHKERRQ(ierr);
560a3df900dSMatthew G Knepley       }
561bdddcaaaSMatthew G Knepley       /* create work vectors for each split */
5622a7a6963SBarry Smith       ierr     = MatCreateVecs(jac->pmat[i],&jac->x[i],&jac->y[i]);CHKERRQ(ierr);
5630298fd71SBarry Smith       ilink->x = jac->x[i]; ilink->y = jac->y[i]; ilink->z = NULL;
564bdddcaaaSMatthew G Knepley       /* compute scatter contexts needed by multiplicative versions and non-default splits */
5650298fd71SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],NULL,&ilink->sctx);CHKERRQ(ierr);
566ed1f0337SMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "nearnullspace", (PetscObject*) &sp);CHKERRQ(ierr);
567ed1f0337SMatthew G Knepley       if (sp) {
568ed1f0337SMatthew G Knepley         ierr = MatSetNearNullSpace(jac->pmat[i], sp);CHKERRQ(ierr);
569ed1f0337SMatthew G Knepley       }
570704ba839SBarry Smith       ilink = ilink->next;
571cf502942SBarry Smith     }
572bdddcaaaSMatthew G Knepley     ierr = VecDestroy(&xtmp);CHKERRQ(ierr);
57397bbdb24SBarry Smith   } else {
574cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
575a3df900dSMatthew G Knepley       Mat pmat;
576a3df900dSMatthew G Knepley 
577a3df900dSMatthew G Knepley       /* Check for preconditioning matrix attached to IS */
578a3df900dSMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &pmat);CHKERRQ(ierr);
579a3df900dSMatthew G Knepley       if (!pmat) {
5805f4ab4e1SJungho Lee         ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
581a3df900dSMatthew G Knepley       }
582704ba839SBarry Smith       ilink = ilink->next;
583cf502942SBarry Smith     }
58497bbdb24SBarry Smith   }
5854e39094bSDmitry Karpeev   if (jac->diag_use_amat) {
586519d70e2SJed Brown     ilink = jac->head;
587519d70e2SJed Brown     if (!jac->mat) {
588785e854fSJed Brown       ierr = PetscMalloc1(nsplit,&jac->mat);CHKERRQ(ierr);
589519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
5905f4ab4e1SJungho Lee         ierr  = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
591519d70e2SJed Brown         ilink = ilink->next;
592519d70e2SJed Brown       }
593519d70e2SJed Brown     } else {
594519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
5955f4ab4e1SJungho Lee         if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);}
596519d70e2SJed Brown         ilink = ilink->next;
597519d70e2SJed Brown       }
598519d70e2SJed Brown     }
599519d70e2SJed Brown   } else {
600519d70e2SJed Brown     jac->mat = jac->pmat;
601519d70e2SJed Brown   }
60297bbdb24SBarry Smith 
60353935eafSBarry Smith   /* Check for null space attached to IS */
60453935eafSBarry Smith   ilink = jac->head;
60553935eafSBarry Smith   for (i=0; i<nsplit; i++) {
60653935eafSBarry Smith     MatNullSpace sp;
60753935eafSBarry Smith 
60853935eafSBarry Smith     ierr = PetscObjectQuery((PetscObject) ilink->is, "nullspace", (PetscObject*) &sp);CHKERRQ(ierr);
60953935eafSBarry Smith     if (sp) {
61053935eafSBarry Smith       ierr = MatSetNullSpace(jac->mat[i], sp);CHKERRQ(ierr);
61153935eafSBarry Smith     }
61253935eafSBarry Smith     ilink = ilink->next;
61353935eafSBarry Smith   }
61453935eafSBarry Smith 
6156c8605c2SJed Brown   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
61668dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
6174e39094bSDmitry Karpeev     /* FIXME: Can/should we reuse jac->mat whenever (jac->diag_use_amat) is true? */
61868dd23aaSBarry Smith     ilink = jac->head;
619e52d2c62SBarry Smith     if (nsplit == 2 && jac->type == PC_COMPOSITE_MULTIPLICATIVE) {
620e52d2c62SBarry Smith       /* 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 */
621e52d2c62SBarry Smith       if (!jac->Afield) {
622e52d2c62SBarry Smith         ierr = PetscCalloc1(nsplit,&jac->Afield);CHKERRQ(ierr);
623e52d2c62SBarry Smith         ierr  = MatGetSubMatrix(pc->mat,ilink->next->is,ilink->is,MAT_INITIAL_MATRIX,&jac->Afield[1]);CHKERRQ(ierr);
624e52d2c62SBarry Smith       } else {
625e52d2c62SBarry Smith         ierr  = MatGetSubMatrix(pc->mat,ilink->next->is,ilink->is,MAT_REUSE_MATRIX,&jac->Afield[1]);CHKERRQ(ierr);
626e52d2c62SBarry Smith       }
627e52d2c62SBarry Smith     } else {
62868dd23aaSBarry Smith       if (!jac->Afield) {
629785e854fSJed Brown         ierr = PetscMalloc1(nsplit,&jac->Afield);CHKERRQ(ierr);
63068dd23aaSBarry Smith         for (i=0; i<nsplit; i++) {
6310298fd71SBarry Smith           ierr  = MatGetSubMatrix(pc->mat,ilink->is,NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
63268dd23aaSBarry Smith           ilink = ilink->next;
63368dd23aaSBarry Smith         }
63468dd23aaSBarry Smith       } else {
63568dd23aaSBarry Smith         for (i=0; i<nsplit; i++) {
6360298fd71SBarry Smith           ierr  = MatGetSubMatrix(pc->mat,ilink->is,NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
63768dd23aaSBarry Smith           ilink = ilink->next;
63868dd23aaSBarry Smith         }
63968dd23aaSBarry Smith       }
64068dd23aaSBarry Smith     }
641e52d2c62SBarry Smith   }
64268dd23aaSBarry Smith 
6433b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
6443b224e63SBarry Smith     IS          ccis;
6454aa3045dSJed Brown     PetscInt    rstart,rend;
646093c86ffSJed Brown     char        lscname[256];
647093c86ffSJed Brown     PetscObject LSC_L;
648ce94432eSBarry Smith 
649ce94432eSBarry Smith     if (nsplit != 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
65068dd23aaSBarry Smith 
651e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
652e6cab6aaSJed Brown     ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
653e6cab6aaSJed Brown 
6543b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
6553b224e63SBarry Smith     if (jac->schur) {
6560298fd71SBarry Smith       KSP kspA = jac->head->ksp, kspInner = NULL, kspUpper = jac->kspupper;
657443836d0SMatthew G Knepley 
658fb3147dbSMatthew G Knepley       ierr  = MatSchurComplementGetKSP(jac->schur, &kspInner);CHKERRQ(ierr);
6593b224e63SBarry Smith       ilink = jac->head;
66049bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
6614e39094bSDmitry Karpeev       if (jac->offdiag_use_amat) {
6624aa3045dSJed Brown 	ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
663c84da90fSDmitry Karpeev       } else {
664df06c9e3SDmitry Karpeev 	ierr  = MatGetSubMatrix(pc->pmat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
665c84da90fSDmitry Karpeev       }
666fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
6673b224e63SBarry Smith       ilink = ilink->next;
66849bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
6694e39094bSDmitry Karpeev       if (jac->offdiag_use_amat) {
6704aa3045dSJed Brown 	ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
671c84da90fSDmitry Karpeev       } else {
672df06c9e3SDmitry Karpeev 	ierr  = MatGetSubMatrix(pc->pmat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
673c84da90fSDmitry Karpeev       }
674fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
675aa6c7ce3SBarry Smith       ierr  = MatSchurComplementUpdateSubMatrices(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);CHKERRQ(ierr);
676a7476a74SDmitry Karpeev       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELFP) {
677a7476a74SDmitry Karpeev 	ierr = MatDestroy(&jac->schurp);CHKERRQ(ierr);
678a7476a74SDmitry Karpeev 	ierr = MatSchurComplementGetPmat(jac->schur,MAT_INITIAL_MATRIX,&jac->schurp);CHKERRQ(ierr);
679a7476a74SDmitry Karpeev       }
680443836d0SMatthew G Knepley       if (kspA != kspInner) {
68123ee1639SBarry Smith         ierr = KSPSetOperators(kspA,jac->mat[0],jac->pmat[0]);CHKERRQ(ierr);
682443836d0SMatthew G Knepley       }
683443836d0SMatthew G Knepley       if (kspUpper != kspA) {
68423ee1639SBarry Smith         ierr = KSPSetOperators(kspUpper,jac->mat[0],jac->pmat[0]);CHKERRQ(ierr);
685443836d0SMatthew G Knepley       }
68623ee1639SBarry Smith       ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac));CHKERRQ(ierr);
6873b224e63SBarry Smith     } else {
688bafc1b83SMatthew G Knepley       const char   *Dprefix;
689470b340bSDmitry Karpeev       char         schurprefix[256], schurmatprefix[256];
690514bf10dSMatthew G Knepley       char         schurtestoption[256];
691bdddcaaaSMatthew G Knepley       MatNullSpace sp;
692514bf10dSMatthew G Knepley       PetscBool    flg;
6933b224e63SBarry Smith 
694a04f6461SBarry Smith       /* extract the A01 and A10 matrices */
6953b224e63SBarry Smith       ilink = jac->head;
69649bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
6974e39094bSDmitry Karpeev       if (jac->offdiag_use_amat) {
6984aa3045dSJed Brown 	ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
699c84da90fSDmitry Karpeev       } else {
700df06c9e3SDmitry Karpeev 	ierr  = MatGetSubMatrix(pc->pmat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
701c84da90fSDmitry Karpeev       }
702fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
7033b224e63SBarry Smith       ilink = ilink->next;
70449bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
7054e39094bSDmitry Karpeev       if (jac->offdiag_use_amat) {
7064aa3045dSJed Brown 	ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
707c84da90fSDmitry Karpeev       } else {
708df06c9e3SDmitry Karpeev 	ierr  = MatGetSubMatrix(pc->pmat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
709c84da90fSDmitry Karpeev       }
710fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
71120252d06SBarry Smith 
712f5236f50SJed Brown       /* Use mat[0] (diagonal block of Amat) preconditioned by pmat[0] to define Schur complement */
71320252d06SBarry Smith       ierr = MatCreate(((PetscObject)jac->mat[0])->comm,&jac->schur);CHKERRQ(ierr);
71420252d06SBarry Smith       ierr = MatSetType(jac->schur,MATSCHURCOMPLEMENT);CHKERRQ(ierr);
715bee83525SDmitry Karpeev       ierr = MatSchurComplementSetSubMatrices(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);CHKERRQ(ierr);
716470b340bSDmitry Karpeev       ierr = PetscSNPrintf(schurmatprefix, sizeof(schurmatprefix), "%sfieldsplit_%s_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr);
717470b340bSDmitry Karpeev       /* Note that the inner KSP is NOT going to inherit this prefix, and if it did, it would be reset just below.  Is that what we want? */
718470b340bSDmitry Karpeev       ierr = MatSetOptionsPrefix(jac->schur,schurmatprefix);CHKERRQ(ierr);
719470b340bSDmitry Karpeev       ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
720895c21f2SBarry Smith       ierr = MatGetNullSpace(jac->mat[1], &sp);CHKERRQ(ierr);
72120252d06SBarry Smith       if (sp) {
72220252d06SBarry Smith         ierr = MatSetNullSpace(jac->schur, sp);CHKERRQ(ierr);
72320252d06SBarry Smith       }
72420252d06SBarry Smith 
72520252d06SBarry Smith       ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_inner_", ilink->splitname);CHKERRQ(ierr);
726c5929fdfSBarry Smith       ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->options,((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);CHKERRQ(ierr);
727514bf10dSMatthew G Knepley       if (flg) {
728514bf10dSMatthew G Knepley         DM  dmInner;
72921635b76SJed Brown         KSP kspInner;
73021635b76SJed Brown 
73121635b76SJed Brown         ierr = MatSchurComplementGetKSP(jac->schur, &kspInner);CHKERRQ(ierr);
73221635b76SJed Brown         ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_inner_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr);
73321635b76SJed Brown         /* Indent this deeper to emphasize the "inner" nature of this solver. */
73421635b76SJed Brown         ierr = PetscObjectIncrementTabLevel((PetscObject)kspInner, (PetscObject) pc, 2);CHKERRQ(ierr);
73521635b76SJed Brown         ierr = KSPSetOptionsPrefix(kspInner, schurprefix);CHKERRQ(ierr);
736514bf10dSMatthew G Knepley 
737514bf10dSMatthew G Knepley         /* Set DM for new solver */
738bafc1b83SMatthew G Knepley         ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr);
73921635b76SJed Brown         ierr = KSPSetDM(kspInner, dmInner);CHKERRQ(ierr);
74021635b76SJed Brown         ierr = KSPSetDMActive(kspInner, PETSC_FALSE);CHKERRQ(ierr);
741514bf10dSMatthew G Knepley       } else {
74221635b76SJed Brown          /* Use the outer solver for the inner solve, but revert the KSPPREONLY from PCFieldSplitSetFields_FieldSplit or
74321635b76SJed Brown           * PCFieldSplitSetIS_FieldSplit. We don't want KSPPREONLY because it makes the Schur complement inexact,
74421635b76SJed Brown           * preventing Schur complement reduction to be an accurate solve. Usually when an iterative solver is used for
74521635b76SJed Brown           * S = D - C A_inner^{-1} B, we expect S to be defined using an accurate definition of A_inner^{-1}, so we make
74621635b76SJed Brown           * GMRES the default. Note that it is also common to use PREONLY for S, in which case S may not be used
74721635b76SJed Brown           * directly, and the user is responsible for setting an inexact method for fieldsplit's A^{-1}. */
74821635b76SJed Brown         ierr = KSPSetType(jac->head->ksp,KSPGMRES);CHKERRQ(ierr);
749514bf10dSMatthew G Knepley         ierr = MatSchurComplementSetKSP(jac->schur,jac->head->ksp);CHKERRQ(ierr);
750bafc1b83SMatthew G Knepley       }
75123ee1639SBarry Smith       ierr = KSPSetOperators(jac->head->ksp,jac->mat[0],jac->pmat[0]);CHKERRQ(ierr);
7525a9f2f41SSatish Balay       ierr = KSPSetFromOptions(jac->head->ksp);CHKERRQ(ierr);
75397bbdb24SBarry Smith       ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
75497bbdb24SBarry Smith 
755443836d0SMatthew G Knepley       ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_upper_", ilink->splitname);CHKERRQ(ierr);
756c5929fdfSBarry Smith       ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->options,((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);CHKERRQ(ierr);
757443836d0SMatthew G Knepley       if (flg) {
758443836d0SMatthew G Knepley         DM dmInner;
759443836d0SMatthew G Knepley 
760443836d0SMatthew G Knepley         ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_upper_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr);
76182f516ccSBarry Smith         ierr = KSPCreate(PetscObjectComm((PetscObject)pc), &jac->kspupper);CHKERRQ(ierr);
762422a814eSBarry Smith         ierr = KSPSetErrorIfNotConverged(jac->kspupper,pc->erroriffailure);CHKERRQ(ierr);
763443836d0SMatthew G Knepley         ierr = KSPSetOptionsPrefix(jac->kspupper, schurprefix);CHKERRQ(ierr);
764443836d0SMatthew G Knepley         ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr);
765443836d0SMatthew G Knepley         ierr = KSPSetDM(jac->kspupper, dmInner);CHKERRQ(ierr);
766443836d0SMatthew G Knepley         ierr = KSPSetDMActive(jac->kspupper, PETSC_FALSE);CHKERRQ(ierr);
767443836d0SMatthew G Knepley         ierr = KSPSetFromOptions(jac->kspupper);CHKERRQ(ierr);
76823ee1639SBarry Smith         ierr = KSPSetOperators(jac->kspupper,jac->mat[0],jac->pmat[0]);CHKERRQ(ierr);
769443836d0SMatthew G Knepley         ierr = VecDuplicate(jac->head->x, &jac->head->z);CHKERRQ(ierr);
770443836d0SMatthew G Knepley       } else {
771443836d0SMatthew G Knepley         jac->kspupper = jac->head->ksp;
772443836d0SMatthew G Knepley         ierr          = PetscObjectReference((PetscObject) jac->head->ksp);CHKERRQ(ierr);
773443836d0SMatthew G Knepley       }
774443836d0SMatthew G Knepley 
775a7476a74SDmitry Karpeev       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELFP) {
776a7476a74SDmitry Karpeev 	ierr = MatSchurComplementGetPmat(jac->schur,MAT_INITIAL_MATRIX,&jac->schurp);CHKERRQ(ierr);
777a7476a74SDmitry Karpeev       }
778ce94432eSBarry Smith       ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&jac->kspschur);CHKERRQ(ierr);
779422a814eSBarry Smith       ierr = KSPSetErrorIfNotConverged(jac->kspschur,pc->erroriffailure);CHKERRQ(ierr);
78097bbdb24SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr);
78120252d06SBarry Smith       ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
78297bbdb24SBarry Smith       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
7837233a360SDmitry Karpeev         PC pcschur;
7847233a360SDmitry Karpeev         ierr = KSPGetPC(jac->kspschur,&pcschur);CHKERRQ(ierr);
7857233a360SDmitry Karpeev         ierr = PCSetType(pcschur,PCNONE);CHKERRQ(ierr);
78697bbdb24SBarry Smith         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
787e74569cdSMatthew G. Knepley       } else if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_FULL) {
788e74569cdSMatthew G. Knepley         ierr = MatSchurComplementComputeExplicitOperator(jac->schur, &jac->schur_user);CHKERRQ(ierr);
78997bbdb24SBarry Smith       }
79023ee1639SBarry Smith       ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac));CHKERRQ(ierr);
79197bbdb24SBarry Smith       ierr = KSPGetOptionsPrefix(jac->head->next->ksp, &Dprefix);CHKERRQ(ierr);
79297bbdb24SBarry Smith       ierr = KSPSetOptionsPrefix(jac->kspschur,         Dprefix);CHKERRQ(ierr);
793b20b4189SMatthew G. Knepley       /* propogate DM */
794b20b4189SMatthew G. Knepley       {
795b20b4189SMatthew G. Knepley         DM sdm;
796b20b4189SMatthew G. Knepley         ierr = KSPGetDM(jac->head->next->ksp, &sdm);CHKERRQ(ierr);
797b20b4189SMatthew G. Knepley         if (sdm) {
798b20b4189SMatthew G. Knepley           ierr = KSPSetDM(jac->kspschur, sdm);CHKERRQ(ierr);
799b20b4189SMatthew G. Knepley           ierr = KSPSetDMActive(jac->kspschur, PETSC_FALSE);CHKERRQ(ierr);
800b20b4189SMatthew G. Knepley         }
801b20b4189SMatthew G. Knepley       }
80297bbdb24SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
80397bbdb24SBarry Smith       /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */
80497bbdb24SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
80597bbdb24SBarry Smith     }
80697bbdb24SBarry Smith 
8075a9f2f41SSatish Balay     /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */
8088caf3d72SBarry Smith     ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_L",ilink->splitname);CHKERRQ(ierr);
809519d70e2SJed Brown     ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
8103b224e63SBarry Smith     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
811c1570756SJed Brown     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);}
8128caf3d72SBarry Smith     ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr);
81397bbdb24SBarry Smith     ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
8145a9f2f41SSatish Balay     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
8150971522cSBarry Smith     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);}
81697bbdb24SBarry Smith   } else {
81768bd789dSDmitry Karpeev     /* set up the individual splits' PCs */
8180971522cSBarry Smith     i     = 0;
8190971522cSBarry Smith     ilink = jac->head;
8200971522cSBarry Smith     while (ilink) {
82123ee1639SBarry Smith       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i]);CHKERRQ(ierr);
8220971522cSBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
8230971522cSBarry Smith       if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);}
8240971522cSBarry Smith       i++;
8250971522cSBarry Smith       ilink = ilink->next;
8260971522cSBarry Smith     }
8273b224e63SBarry Smith   }
8283b224e63SBarry Smith 
829c1570756SJed Brown   jac->suboptionsset = PETSC_TRUE;
8300971522cSBarry Smith   PetscFunctionReturn(0);
8310971522cSBarry Smith }
8320971522cSBarry Smith 
8335a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
834ca9f406cSSatish Balay   (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
835ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
8369df09d43SBarry Smith    PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL) ||\
8375a9f2f41SSatish Balay    KSPSolve(ilink->ksp,ilink->x,ilink->y) ||                               \
8389df09d43SBarry Smith    PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL) ||\
839ca9f406cSSatish Balay    VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) ||  \
840ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
84179416396SBarry Smith 
8420971522cSBarry Smith #undef __FUNCT__
8433b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur"
8443b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
8453b224e63SBarry Smith {
8463b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
8473b224e63SBarry Smith   PetscErrorCode    ierr;
8483b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
849443836d0SMatthew G Knepley   KSP               kspA   = ilinkA->ksp, kspLower = kspA, kspUpper = jac->kspupper;
8503b224e63SBarry Smith 
8513b224e63SBarry Smith   PetscFunctionBegin;
852c5d2311dSJed Brown   switch (jac->schurfactorization) {
853c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_DIAG:
854a04f6461SBarry Smith     /* [A00 0; 0 -S], positive definite, suitable for MINRES */
855c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
856c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
857c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8589df09d43SBarry Smith     ierr = PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
859443836d0SMatthew G Knepley     ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
8609df09d43SBarry Smith     ierr = PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
861c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
862c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8639df09d43SBarry Smith     ierr = PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
864c5d2311dSJed Brown     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
8659df09d43SBarry Smith     ierr = PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
866c5d2311dSJed Brown     ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr);
867c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
868c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
869c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
870c5d2311dSJed Brown     break;
871c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_LOWER:
872a04f6461SBarry Smith     /* [A00 0; A10 S], suitable for left preconditioning */
873c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
874c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8759df09d43SBarry Smith     ierr = PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
876443836d0SMatthew G Knepley     ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
8779df09d43SBarry Smith     ierr = PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
878c5d2311dSJed Brown     ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
879c5d2311dSJed Brown     ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr);
880c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
881c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
882c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8839df09d43SBarry Smith     ierr = PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
884c5d2311dSJed Brown     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
8859df09d43SBarry Smith     ierr = PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
886c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
887c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
888c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
889c5d2311dSJed Brown     break;
890c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_UPPER:
891a04f6461SBarry Smith     /* [A00 A01; 0 S], suitable for right preconditioning */
892c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
893c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8949df09d43SBarry Smith     ierr = PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
895c5d2311dSJed Brown     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
8969df09d43SBarry Smith     ierr = PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);    ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
897c5d2311dSJed Brown     ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr);
898c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
899c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
900c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9019df09d43SBarry Smith     ierr = PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
902443836d0SMatthew G Knepley     ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
9039df09d43SBarry Smith     ierr = PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
904c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
905c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
906c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
907c5d2311dSJed Brown     break;
908c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_FULL:
909a04f6461SBarry Smith     /* [1 0; A10 A00^{-1} 1] [A00 0; 0 S] [1 A00^{-1}A01; 0 1], an exact solve if applied exactly, needs one extra solve with A */
9103b224e63SBarry Smith     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9113b224e63SBarry Smith     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9129df09d43SBarry Smith     ierr = PetscLogEventBegin(KSP_Solve_FS_L,kspLower,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
913443836d0SMatthew G Knepley     ierr = KSPSolve(kspLower,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
9149df09d43SBarry Smith     ierr = PetscLogEventEnd(KSP_Solve_FS_L,kspLower,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
9153b224e63SBarry Smith     ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
9163b224e63SBarry Smith     ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
9173b224e63SBarry Smith     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9183b224e63SBarry Smith     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9193b224e63SBarry Smith 
9209df09d43SBarry Smith     ierr = PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
9213b224e63SBarry Smith     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
9229df09d43SBarry Smith     ierr = PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
9233b224e63SBarry Smith     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9243b224e63SBarry Smith     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9253b224e63SBarry Smith 
926443836d0SMatthew G Knepley     if (kspUpper == kspA) {
9273b224e63SBarry Smith       ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
9283b224e63SBarry Smith       ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
9299df09d43SBarry Smith       ierr = PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
930443836d0SMatthew G Knepley       ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
9319df09d43SBarry Smith       ierr = PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
932443836d0SMatthew G Knepley     } else {
9339df09d43SBarry Smith       ierr = PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
934443836d0SMatthew G Knepley       ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
9359df09d43SBarry Smith       ierr = PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
936443836d0SMatthew G Knepley       ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
9379df09d43SBarry Smith       ierr = PetscLogEventBegin(KSP_Solve_FS_U,kspUpper,ilinkA->x,ilinkA->z,NULL);CHKERRQ(ierr);
938443836d0SMatthew G Knepley       ierr = KSPSolve(kspUpper,ilinkA->x,ilinkA->z);CHKERRQ(ierr);
9399df09d43SBarry Smith       ierr = PetscLogEventEnd(KSP_Solve_FS_U,kspUpper,ilinkA->x,ilinkA->z,NULL);CHKERRQ(ierr);
940443836d0SMatthew G Knepley       ierr = VecAXPY(ilinkA->y,-1.0,ilinkA->z);CHKERRQ(ierr);
941443836d0SMatthew G Knepley     }
9423b224e63SBarry Smith     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9433b224e63SBarry Smith     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
944c5d2311dSJed Brown   }
9453b224e63SBarry Smith   PetscFunctionReturn(0);
9463b224e63SBarry Smith }
9473b224e63SBarry Smith 
9483b224e63SBarry Smith #undef __FUNCT__
9490971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
9500971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
9510971522cSBarry Smith {
9520971522cSBarry Smith   PC_FieldSplit      *jac = (PC_FieldSplit*)pc->data;
9530971522cSBarry Smith   PetscErrorCode     ierr;
9545a9f2f41SSatish Balay   PC_FieldSplitLink  ilink = jac->head;
955939b8a20SBarry Smith   PetscInt           cnt,bs;
956568ed31eSHong Zhang   KSPConvergedReason reason;
9570971522cSBarry Smith 
9580971522cSBarry Smith   PetscFunctionBegin;
95979416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
9601b9fc7fcSBarry Smith     if (jac->defaultsplit) {
961939b8a20SBarry Smith       ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
962ce94432eSBarry Smith       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);
963939b8a20SBarry Smith       ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr);
964ce94432eSBarry Smith       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);
9650971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
9665a9f2f41SSatish Balay       while (ilink) {
9679df09d43SBarry Smith         ierr  = PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
9685a9f2f41SSatish Balay         ierr  = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
9699df09d43SBarry Smith         ierr  = PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
970568ed31eSHong Zhang         ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
971568ed31eSHong Zhang         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
972568ed31eSHong Zhang           pc->failedreason = PC_SUBPC_ERROR;
973568ed31eSHong Zhang         }
9745a9f2f41SSatish Balay         ilink = ilink->next;
9750971522cSBarry Smith       }
9760971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
9771b9fc7fcSBarry Smith     } else {
978efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
9795a9f2f41SSatish Balay       while (ilink) {
9805a9f2f41SSatish Balay         ierr  = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
981568ed31eSHong Zhang         ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
982568ed31eSHong Zhang         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
983568ed31eSHong Zhang           pc->failedreason = PC_SUBPC_ERROR;
984568ed31eSHong Zhang         }
9855a9f2f41SSatish Balay         ilink = ilink->next;
9861b9fc7fcSBarry Smith       }
9871b9fc7fcSBarry Smith     }
988e52d2c62SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE && jac->nsplits == 2) {
989e52d2c62SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
990e52d2c62SBarry Smith     /* solve on first block for first block variables */
991e52d2c62SBarry Smith     ierr = VecScatterBegin(ilink->sctx,x,ilink->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
992e52d2c62SBarry Smith     ierr = VecScatterEnd(ilink->sctx,x,ilink->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9939df09d43SBarry Smith     ierr = PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
994e52d2c62SBarry Smith     ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
9959df09d43SBarry Smith     ierr = PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
996568ed31eSHong Zhang     ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
997568ed31eSHong Zhang     if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
998568ed31eSHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
999568ed31eSHong Zhang     }
1000e52d2c62SBarry Smith     ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1001e52d2c62SBarry Smith     ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1002e52d2c62SBarry Smith 
1003e52d2c62SBarry Smith     /* compute the residual only onto second block variables using first block variables */
1004e52d2c62SBarry Smith     ierr = MatMult(jac->Afield[1],ilink->y,ilink->next->x);CHKERRQ(ierr);
1005e52d2c62SBarry Smith     ilink = ilink->next;
1006e52d2c62SBarry Smith     ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
1007e52d2c62SBarry Smith     ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1008e52d2c62SBarry Smith     ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1009e52d2c62SBarry Smith 
1010e52d2c62SBarry Smith     /* solve on second block variables */
10119df09d43SBarry Smith     ierr = PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
1012e52d2c62SBarry Smith     ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
10139df09d43SBarry Smith     ierr = PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
1014568ed31eSHong Zhang     ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
1015568ed31eSHong Zhang     if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1016568ed31eSHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
1017568ed31eSHong Zhang     }
1018e52d2c62SBarry Smith     ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1019e52d2c62SBarry Smith     ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
102016913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
102179416396SBarry Smith     if (!jac->w1) {
102279416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
102379416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
102479416396SBarry Smith     }
1025efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
10265a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
1027568ed31eSHong Zhang     ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
1028568ed31eSHong Zhang     if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1029568ed31eSHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
1030568ed31eSHong Zhang     }
10313e197d65SBarry Smith     cnt  = 1;
10325a9f2f41SSatish Balay     while (ilink->next) {
10335a9f2f41SSatish Balay       ilink = ilink->next;
10343e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
10353e197d65SBarry Smith       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
10363e197d65SBarry Smith       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
10373e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10383e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10399df09d43SBarry Smith       ierr = PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
10403e197d65SBarry Smith       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
10419df09d43SBarry Smith       ierr = PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
10422618eb8fSHong Zhang       ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
10432618eb8fSHong Zhang       if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
10442618eb8fSHong Zhang         pc->failedreason = PC_SUBPC_ERROR;
10452618eb8fSHong Zhang       }
10463e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10473e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10483e197d65SBarry Smith     }
104951f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
105011755939SBarry Smith       cnt -= 2;
105151f519a2SBarry Smith       while (ilink->previous) {
105251f519a2SBarry Smith         ilink = ilink->previous;
105311755939SBarry Smith         /* compute the residual only over the part of the vector needed */
105411755939SBarry Smith         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
105511755939SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
105611755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
105711755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10589df09d43SBarry Smith         ierr = PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
105911755939SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
10609df09d43SBarry Smith         ierr = PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
1061568ed31eSHong Zhang         ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
1062568ed31eSHong Zhang         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1063568ed31eSHong Zhang           pc->failedreason = PC_SUBPC_ERROR;
1064568ed31eSHong Zhang         }
106511755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
106611755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
106751f519a2SBarry Smith       }
106811755939SBarry Smith     }
1069ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
10700971522cSBarry Smith   PetscFunctionReturn(0);
10710971522cSBarry Smith }
10720971522cSBarry Smith 
1073421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
1074ca9f406cSSatish Balay   (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
1075ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
10769df09d43SBarry Smith    PetscLogEventBegin(ilink->event,ilink->ksp,ilink->y,ilink->x,NULL) || \
1077421e10b8SBarry Smith    KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) ||                  \
10789df09d43SBarry Smith    PetscLogEventBegin(ilink->event,ilink->ksp,ilink->y,ilink->x,NULL) || \
1079ca9f406cSSatish Balay    VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
1080ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
1081421e10b8SBarry Smith 
1082421e10b8SBarry Smith #undef __FUNCT__
10838c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit"
1084421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
1085421e10b8SBarry Smith {
1086421e10b8SBarry Smith   PC_FieldSplit      *jac = (PC_FieldSplit*)pc->data;
1087421e10b8SBarry Smith   PetscErrorCode     ierr;
1088421e10b8SBarry Smith   PC_FieldSplitLink  ilink = jac->head;
1089939b8a20SBarry Smith   PetscInt           bs;
1090291dd214SHong Zhang   KSPConvergedReason reason;
1091421e10b8SBarry Smith 
1092421e10b8SBarry Smith   PetscFunctionBegin;
1093421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
1094421e10b8SBarry Smith     if (jac->defaultsplit) {
1095939b8a20SBarry Smith       ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
1096ce94432eSBarry Smith       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);
1097939b8a20SBarry Smith       ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr);
1098ce94432eSBarry Smith       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);
1099421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
1100421e10b8SBarry Smith       while (ilink) {
11019df09d43SBarry Smith         ierr = PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
1102421e10b8SBarry Smith         ierr  = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
11039df09d43SBarry Smith         ierr = PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
11042618eb8fSHong Zhang         ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
11052618eb8fSHong Zhang         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
11062618eb8fSHong Zhang           pc->failedreason = PC_SUBPC_ERROR;
11072618eb8fSHong Zhang         }
1108421e10b8SBarry Smith         ilink = ilink->next;
1109421e10b8SBarry Smith       }
1110421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
1111421e10b8SBarry Smith     } else {
1112421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
1113421e10b8SBarry Smith       while (ilink) {
1114421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
11152618eb8fSHong Zhang         ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
11162618eb8fSHong Zhang         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
11172618eb8fSHong Zhang           pc->failedreason = PC_SUBPC_ERROR;
11182618eb8fSHong Zhang         }
1119421e10b8SBarry Smith         ilink = ilink->next;
1120421e10b8SBarry Smith       }
1121421e10b8SBarry Smith     }
1122421e10b8SBarry Smith   } else {
1123421e10b8SBarry Smith     if (!jac->w1) {
1124421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
1125421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
1126421e10b8SBarry Smith     }
1127421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
1128421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
1129421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
11302618eb8fSHong Zhang       ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
11312618eb8fSHong Zhang       if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
11322618eb8fSHong Zhang         pc->failedreason = PC_SUBPC_ERROR;
11332618eb8fSHong Zhang       }
1134421e10b8SBarry Smith       while (ilink->next) {
1135421e10b8SBarry Smith         ilink = ilink->next;
11369989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
1137421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
1138421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
1139421e10b8SBarry Smith       }
1140421e10b8SBarry Smith       while (ilink->previous) {
1141421e10b8SBarry Smith         ilink = ilink->previous;
11429989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
1143421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
1144421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
1145421e10b8SBarry Smith       }
1146421e10b8SBarry Smith     } else {
1147421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
1148421e10b8SBarry Smith         ilink = ilink->next;
1149421e10b8SBarry Smith       }
1150421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
11512618eb8fSHong Zhang       ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
11522618eb8fSHong Zhang       if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
11532618eb8fSHong Zhang         pc->failedreason = PC_SUBPC_ERROR;
11542618eb8fSHong Zhang       }
1155421e10b8SBarry Smith       while (ilink->previous) {
1156421e10b8SBarry Smith         ilink = ilink->previous;
11579989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
1158421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
1159421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
1160421e10b8SBarry Smith       }
1161421e10b8SBarry Smith     }
1162421e10b8SBarry Smith   }
1163421e10b8SBarry Smith   PetscFunctionReturn(0);
1164421e10b8SBarry Smith }
1165421e10b8SBarry Smith 
11660971522cSBarry Smith #undef __FUNCT__
1167574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit"
1168574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc)
11690971522cSBarry Smith {
11700971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
11710971522cSBarry Smith   PetscErrorCode    ierr;
11725a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
11730971522cSBarry Smith 
11740971522cSBarry Smith   PetscFunctionBegin;
11755a9f2f41SSatish Balay   while (ilink) {
1176574deadeSBarry Smith     ierr  = KSPReset(ilink->ksp);CHKERRQ(ierr);
1177fcfd50ebSBarry Smith     ierr  = VecDestroy(&ilink->x);CHKERRQ(ierr);
1178fcfd50ebSBarry Smith     ierr  = VecDestroy(&ilink->y);CHKERRQ(ierr);
1179443836d0SMatthew G Knepley     ierr  = VecDestroy(&ilink->z);CHKERRQ(ierr);
1180fcfd50ebSBarry Smith     ierr  = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr);
11811c7cfba8SBarry Smith     if (!ilink->is_orig) {             /* save the original IS */
11826dbb499eSCian Wilson       ierr = PetscObjectReference((PetscObject)ilink->is);CHKERRQ(ierr);
11836dbb499eSCian Wilson       ilink->is_orig = ilink->is;
11846dbb499eSCian Wilson     }
1185fcfd50ebSBarry Smith     ierr  = ISDestroy(&ilink->is);CHKERRQ(ierr);
1186b5787286SJed Brown     ierr  = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
11875a9f2f41SSatish Balay     next  = ilink->next;
11885a9f2f41SSatish Balay     ilink = next;
11890971522cSBarry Smith   }
119005b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
1191574deadeSBarry Smith   if (jac->mat && jac->mat != jac->pmat) {
1192574deadeSBarry Smith     ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);
1193574deadeSBarry Smith   } else if (jac->mat) {
11940298fd71SBarry Smith     jac->mat = NULL;
1195574deadeSBarry Smith   }
119697bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
119768dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
11986bf464f9SBarry Smith   ierr       = VecDestroy(&jac->w1);CHKERRQ(ierr);
11996bf464f9SBarry Smith   ierr       = VecDestroy(&jac->w2);CHKERRQ(ierr);
12006bf464f9SBarry Smith   ierr       = MatDestroy(&jac->schur);CHKERRQ(ierr);
1201a7476a74SDmitry Karpeev   ierr       = MatDestroy(&jac->schurp);CHKERRQ(ierr);
12026bf464f9SBarry Smith   ierr       = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
12036bf464f9SBarry Smith   ierr       = KSPDestroy(&jac->kspschur);CHKERRQ(ierr);
1204d78dad28SBarry Smith   ierr       = KSPDestroy(&jac->kspupper);CHKERRQ(ierr);
12056bf464f9SBarry Smith   ierr       = MatDestroy(&jac->B);CHKERRQ(ierr);
12066bf464f9SBarry Smith   ierr       = MatDestroy(&jac->C);CHKERRQ(ierr);
120763ec74ffSBarry Smith   jac->reset = PETSC_TRUE;
12086dbb499eSCian Wilson   jac->isrestrict = PETSC_FALSE;
1209574deadeSBarry Smith   PetscFunctionReturn(0);
1210574deadeSBarry Smith }
1211574deadeSBarry Smith 
1212574deadeSBarry Smith #undef __FUNCT__
1213574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
1214574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
1215574deadeSBarry Smith {
1216574deadeSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1217574deadeSBarry Smith   PetscErrorCode    ierr;
1218574deadeSBarry Smith   PC_FieldSplitLink ilink = jac->head,next;
1219574deadeSBarry Smith 
1220574deadeSBarry Smith   PetscFunctionBegin;
1221574deadeSBarry Smith   ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr);
1222574deadeSBarry Smith   while (ilink) {
12236bf464f9SBarry Smith     ierr  = KSPDestroy(&ilink->ksp);CHKERRQ(ierr);
12246dbb499eSCian Wilson     ierr  = ISDestroy(&ilink->is_orig);CHKERRQ(ierr);
1225574deadeSBarry Smith     next  = ilink->next;
1226574deadeSBarry Smith     ierr  = PetscFree(ilink->splitname);CHKERRQ(ierr);
1227574deadeSBarry Smith     ierr  = PetscFree(ilink->fields);CHKERRQ(ierr);
12285d4c12cdSJungho Lee     ierr  = PetscFree(ilink->fields_col);CHKERRQ(ierr);
1229574deadeSBarry Smith     ierr  = PetscFree(ilink);CHKERRQ(ierr);
1230574deadeSBarry Smith     ilink = next;
1231574deadeSBarry Smith   }
1232574deadeSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
1233c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1234bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",NULL);CHKERRQ(ierr);
1235bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",NULL);CHKERRQ(ierr);
1236bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",NULL);CHKERRQ(ierr);
1237bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",NULL);CHKERRQ(ierr);
1238bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",NULL);CHKERRQ(ierr);
123929f8a81cSJed Brown   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",NULL);CHKERRQ(ierr);
124037a82bf0SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",NULL);CHKERRQ(ierr);
1241bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",NULL);CHKERRQ(ierr);
12426dbb499eSCian Wilson   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitRestrictIS_C",NULL);CHKERRQ(ierr);
12430971522cSBarry Smith   PetscFunctionReturn(0);
12440971522cSBarry Smith }
12450971522cSBarry Smith 
12460971522cSBarry Smith #undef __FUNCT__
12470971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
12484416b707SBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PetscOptionItems *PetscOptionsObject,PC pc)
12490971522cSBarry Smith {
12501b9fc7fcSBarry Smith   PetscErrorCode  ierr;
12516c924f48SJed Brown   PetscInt        bs;
1252bc59fbc5SBarry Smith   PetscBool       flg,stokes = PETSC_FALSE;
12539dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
12543b224e63SBarry Smith   PCCompositeType ctype;
12551b9fc7fcSBarry Smith 
12560971522cSBarry Smith   PetscFunctionBegin;
1257e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"FieldSplit options");CHKERRQ(ierr);
12584ab8060aSDmitry Karpeev   ierr = PetscOptionsBool("-pc_fieldsplit_dm_splits","Whether to use DMCreateFieldDecomposition() for splits","PCFieldSplitSetDMSplits",jac->dm_splits,&jac->dm_splits,NULL);CHKERRQ(ierr);
125951f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
126051f519a2SBarry Smith   if (flg) {
126151f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
126251f519a2SBarry Smith   }
12632686e3e9SMatthew G. Knepley   jac->diag_use_amat = pc->useAmat;
12642686e3e9SMatthew G. Knepley   ierr = 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);CHKERRQ(ierr);
12652686e3e9SMatthew G. Knepley   jac->offdiag_use_amat = pc->useAmat;
12662686e3e9SMatthew G. Knepley   ierr = 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);CHKERRQ(ierr);
126700216929SDmitry Karpeev   /* FIXME: No programmatic equivalent to the following. */
1268c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);CHKERRQ(ierr);
1269c0adfefeSBarry Smith   if (stokes) {
1270c0adfefeSBarry Smith     ierr          = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
1271c0adfefeSBarry Smith     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
1272c0adfefeSBarry Smith   }
1273c0adfefeSBarry Smith 
12743b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
12753b224e63SBarry Smith   if (flg) {
12763b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
12773b224e63SBarry Smith   }
1278c30613efSMatthew Knepley   /* Only setup fields once */
1279c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
1280d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
1281d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
12826c924f48SJed Brown     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
12836c924f48SJed Brown     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
1284d32f9abdSBarry Smith   }
1285c5d2311dSJed Brown   if (jac->type == PC_COMPOSITE_SCHUR) {
1286c5929fdfSBarry Smith     ierr = PetscOptionsGetEnum(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr);
1287c9c6ffaaSJed Brown     if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);}
12880298fd71SBarry Smith     ierr = 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);CHKERRQ(ierr);
128929f8a81cSJed Brown     ierr = PetscOptionsEnum("-pc_fieldsplit_schur_precondition","How to build preconditioner for Schur complement","PCFieldSplitSetSchurPre",PCFieldSplitSchurPreTypes,(PetscEnum)jac->schurpre,(PetscEnum*)&jac->schurpre,NULL);CHKERRQ(ierr);
1290c5d2311dSJed Brown   }
12911b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
12920971522cSBarry Smith   PetscFunctionReturn(0);
12930971522cSBarry Smith }
12940971522cSBarry Smith 
12950971522cSBarry Smith /*------------------------------------------------------------------------------------*/
12960971522cSBarry Smith 
12970971522cSBarry Smith #undef __FUNCT__
12980971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
12991e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
13000971522cSBarry Smith {
130197bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
13020971522cSBarry Smith   PetscErrorCode    ierr;
13035a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
130469a612a9SBarry Smith   char              prefix[128];
13055d4c12cdSJungho Lee   PetscInt          i;
13060971522cSBarry Smith 
13070971522cSBarry Smith   PetscFunctionBegin;
13086c924f48SJed Brown   if (jac->splitdefined) {
13096c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
13106c924f48SJed Brown     PetscFunctionReturn(0);
13116c924f48SJed Brown   }
131251f519a2SBarry Smith   for (i=0; i<n; i++) {
1313e32f2f54SBarry Smith     if (fields[i] >= jac->bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs);
1314e32f2f54SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
131551f519a2SBarry Smith   }
1316b00a9115SJed Brown   ierr = PetscNew(&ilink);CHKERRQ(ierr);
1317a04f6461SBarry Smith   if (splitname) {
1318db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1319a04f6461SBarry Smith   } else {
1320785e854fSJed Brown     ierr = PetscMalloc1(3,&ilink->splitname);CHKERRQ(ierr);
1321a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
1322a04f6461SBarry Smith   }
13239df09d43SBarry Smith   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 */
1324785e854fSJed Brown   ierr = PetscMalloc1(n,&ilink->fields);CHKERRQ(ierr);
13255a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
1326785e854fSJed Brown   ierr = PetscMalloc1(n,&ilink->fields_col);CHKERRQ(ierr);
13275d4c12cdSJungho Lee   ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr);
13282fa5cd67SKarl Rupp 
13295a9f2f41SSatish Balay   ilink->nfields = n;
13300298fd71SBarry Smith   ilink->next    = NULL;
1331ce94432eSBarry Smith   ierr           = KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);CHKERRQ(ierr);
1332422a814eSBarry Smith   ierr           = KSPSetErrorIfNotConverged(ilink->ksp,pc->erroriffailure);CHKERRQ(ierr);
133320252d06SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
13345a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
13359005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
133669a612a9SBarry Smith 
13378caf3d72SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);CHKERRQ(ierr);
13385a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
13390971522cSBarry Smith 
13400971522cSBarry Smith   if (!next) {
13415a9f2f41SSatish Balay     jac->head       = ilink;
13420298fd71SBarry Smith     ilink->previous = NULL;
13430971522cSBarry Smith   } else {
13440971522cSBarry Smith     while (next->next) {
13450971522cSBarry Smith       next = next->next;
13460971522cSBarry Smith     }
13475a9f2f41SSatish Balay     next->next      = ilink;
134851f519a2SBarry Smith     ilink->previous = next;
13490971522cSBarry Smith   }
13500971522cSBarry Smith   jac->nsplits++;
13510971522cSBarry Smith   PetscFunctionReturn(0);
13520971522cSBarry Smith }
13530971522cSBarry Smith 
1354e69d4d44SBarry Smith #undef __FUNCT__
1355e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
13561e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
1357e69d4d44SBarry Smith {
1358e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1359e69d4d44SBarry Smith   PetscErrorCode ierr;
1360e69d4d44SBarry Smith 
1361e69d4d44SBarry Smith   PetscFunctionBegin;
1362785e854fSJed Brown   ierr = PetscMalloc1(jac->nsplits,subksp);CHKERRQ(ierr);
1363e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
13642fa5cd67SKarl Rupp 
1365e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
136613e0d083SBarry Smith   if (n) *n = jac->nsplits;
1367e69d4d44SBarry Smith   PetscFunctionReturn(0);
1368e69d4d44SBarry Smith }
13690971522cSBarry Smith 
13700971522cSBarry Smith #undef __FUNCT__
137169a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
13721e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
13730971522cSBarry Smith {
13740971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
13750971522cSBarry Smith   PetscErrorCode    ierr;
13760971522cSBarry Smith   PetscInt          cnt   = 0;
13775a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
13780971522cSBarry Smith 
13790971522cSBarry Smith   PetscFunctionBegin;
1380785e854fSJed Brown   ierr = PetscMalloc1(jac->nsplits,subksp);CHKERRQ(ierr);
13815a9f2f41SSatish Balay   while (ilink) {
13825a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
13835a9f2f41SSatish Balay     ilink            = ilink->next;
13840971522cSBarry Smith   }
13855d480477SMatthew G Knepley   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);
138613e0d083SBarry Smith   if (n) *n = jac->nsplits;
13870971522cSBarry Smith   PetscFunctionReturn(0);
13880971522cSBarry Smith }
13890971522cSBarry Smith 
1390704ba839SBarry Smith #undef __FUNCT__
13916dbb499eSCian Wilson #define __FUNCT__ "PCFieldSplitRestrictIS"
13926dbb499eSCian Wilson /*@C
13936dbb499eSCian Wilson     PCFieldSplitRestrictIS - Restricts the fieldsplit ISs to be within a given IS.
13946dbb499eSCian Wilson 
13956dbb499eSCian Wilson     Input Parameters:
13966dbb499eSCian Wilson +   pc  - the preconditioner context
13976dbb499eSCian Wilson +   is - the index set that defines the indices to which the fieldsplit is to be restricted
13986dbb499eSCian Wilson 
13996dbb499eSCian Wilson     Level: advanced
14006dbb499eSCian Wilson 
14016dbb499eSCian Wilson @*/
14026dbb499eSCian Wilson PetscErrorCode  PCFieldSplitRestrictIS(PC pc,IS isy)
14036dbb499eSCian Wilson {
14046dbb499eSCian Wilson   PetscErrorCode ierr;
14056dbb499eSCian Wilson 
14066dbb499eSCian Wilson   PetscFunctionBegin;
14076dbb499eSCian Wilson   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
14086dbb499eSCian Wilson   PetscValidHeaderSpecific(isy,IS_CLASSID,2);
14090246f55cSBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitRestrictIS_C",(PC,IS),(pc,isy));CHKERRQ(ierr);
14106dbb499eSCian Wilson   PetscFunctionReturn(0);
14116dbb499eSCian Wilson }
14126dbb499eSCian Wilson 
14136dbb499eSCian Wilson 
14146dbb499eSCian Wilson #undef __FUNCT__
14156dbb499eSCian Wilson #define __FUNCT__ "PCFieldSplitRestrictIS_FieldSplit"
14166dbb499eSCian Wilson static PetscErrorCode  PCFieldSplitRestrictIS_FieldSplit(PC pc, IS isy)
14176dbb499eSCian Wilson {
14186dbb499eSCian Wilson   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
14196dbb499eSCian Wilson   PetscErrorCode    ierr;
14206dbb499eSCian Wilson   PC_FieldSplitLink ilink = jac->head, next;
14216dbb499eSCian Wilson   PetscInt          localsize,size,sizez,i;
14226dbb499eSCian Wilson   const PetscInt    *ind, *indz;
14236dbb499eSCian Wilson   PetscInt          *indc, *indcz;
14246dbb499eSCian Wilson   PetscBool         flg;
14256dbb499eSCian Wilson 
14266dbb499eSCian Wilson   PetscFunctionBegin;
14276dbb499eSCian Wilson   ierr = ISGetLocalSize(isy,&localsize);CHKERRQ(ierr);
14286dbb499eSCian Wilson   ierr = MPI_Scan(&localsize,&size,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)isy));CHKERRQ(ierr);
14296dbb499eSCian Wilson   size -= localsize;
14306dbb499eSCian Wilson   while(ilink) {
14316dbb499eSCian Wilson     IS isrl,isr;
14321c7cfba8SBarry Smith     PC subpc;
14336dbb499eSCian Wilson     if (jac->reset) {
14346dbb499eSCian Wilson       ierr = ISEmbed(ilink->is_orig, isy, PETSC_TRUE, &isrl);CHKERRQ(ierr);
14356dbb499eSCian Wilson     } else {
14366dbb499eSCian Wilson       ierr = ISEmbed(ilink->is, isy, PETSC_TRUE, &isrl);CHKERRQ(ierr);
14376dbb499eSCian Wilson     }
14386dbb499eSCian Wilson     ierr          = ISGetLocalSize(isrl,&localsize);CHKERRQ(ierr);
14396dbb499eSCian Wilson     ierr          = PetscMalloc1(localsize,&indc);CHKERRQ(ierr);
14406dbb499eSCian Wilson     ierr          = ISGetIndices(isrl,&ind);CHKERRQ(ierr);
14416dbb499eSCian Wilson     ierr          = PetscMemcpy(indc,ind,localsize*sizeof(PetscInt));CHKERRQ(ierr);
14426dbb499eSCian Wilson     ierr          = ISRestoreIndices(isrl,&ind);CHKERRQ(ierr);
14436dbb499eSCian Wilson     ierr          = ISDestroy(&isrl);CHKERRQ(ierr);
14446dbb499eSCian Wilson     for (i=0; i<localsize; i++) *(indc+i) += size;
14456dbb499eSCian Wilson     ierr          = ISCreateGeneral(PetscObjectComm((PetscObject)isy),localsize,indc,PETSC_OWN_POINTER,&isr);CHKERRQ(ierr);
14466dbb499eSCian Wilson     ierr          = PetscObjectReference((PetscObject)isr);CHKERRQ(ierr);
14476dbb499eSCian Wilson     ierr          = ISDestroy(&ilink->is);CHKERRQ(ierr);
14486dbb499eSCian Wilson     ilink->is     = isr;
14496dbb499eSCian Wilson     ierr          = PetscObjectReference((PetscObject)isr);CHKERRQ(ierr);
14506dbb499eSCian Wilson     ierr          = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
14516dbb499eSCian Wilson     ilink->is_col = isr;
14526dbb499eSCian Wilson     ierr          = ISDestroy(&isr);CHKERRQ(ierr);
14536dbb499eSCian Wilson     ierr          = KSPGetPC(ilink->ksp, &subpc);CHKERRQ(ierr);
14546dbb499eSCian Wilson     ierr          = PetscObjectTypeCompare((PetscObject)subpc,PCFIELDSPLIT,&flg);CHKERRQ(ierr);
14556dbb499eSCian Wilson     if(flg) {
14566dbb499eSCian Wilson       IS iszl,isz;
14576dbb499eSCian Wilson       MPI_Comm comm;
14586dbb499eSCian Wilson       if (jac->reset) {
14596dbb499eSCian Wilson         ierr = ISGetLocalSize(ilink->is_orig,&localsize);CHKERRQ(ierr);
14606dbb499eSCian Wilson         comm = PetscObjectComm((PetscObject)ilink->is_orig);
14616dbb499eSCian Wilson         ierr = ISEmbed(isy, ilink->is_orig, PETSC_TRUE, &iszl);CHKERRQ(ierr);
14626dbb499eSCian Wilson       } else {
14636dbb499eSCian Wilson         ierr = ISGetLocalSize(ilink->is,&localsize);CHKERRQ(ierr);
14646dbb499eSCian Wilson         comm = PetscObjectComm((PetscObject)ilink->is);
14656dbb499eSCian Wilson         ierr = ISEmbed(isy, ilink->is, PETSC_TRUE, &iszl);CHKERRQ(ierr);
14666dbb499eSCian Wilson       }
14676dbb499eSCian Wilson       ierr   = MPI_Scan(&localsize,&sizez,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr);
14686dbb499eSCian Wilson       sizez -= localsize;
14696dbb499eSCian Wilson       ierr   = ISGetLocalSize(iszl,&localsize);CHKERRQ(ierr);
14706dbb499eSCian Wilson       ierr   = PetscMalloc1(localsize,&indcz);CHKERRQ(ierr);
14716dbb499eSCian Wilson       ierr   = ISGetIndices(iszl,&indz);CHKERRQ(ierr);
14726dbb499eSCian Wilson       ierr   = PetscMemcpy(indcz,indz,localsize*sizeof(PetscInt));CHKERRQ(ierr);
14736dbb499eSCian Wilson       ierr   = ISRestoreIndices(iszl,&indz);CHKERRQ(ierr);
14746dbb499eSCian Wilson       ierr   = ISDestroy(&iszl);CHKERRQ(ierr);
14756dbb499eSCian Wilson       for (i=0; i<localsize; i++) *(indcz+i) += sizez;
14766dbb499eSCian Wilson       ierr   = ISCreateGeneral(comm,localsize,indcz,PETSC_OWN_POINTER,&isz);CHKERRQ(ierr);
14776dbb499eSCian Wilson       ierr   = PCFieldSplitRestrictIS(subpc,isz);CHKERRQ(ierr);
14786dbb499eSCian Wilson       ierr   = ISDestroy(&isz);CHKERRQ(ierr);
14796dbb499eSCian Wilson     }
14806dbb499eSCian Wilson     next = ilink->next;
14816dbb499eSCian Wilson     ilink = next;
14826dbb499eSCian Wilson   }
14836dbb499eSCian Wilson   jac->isrestrict = PETSC_TRUE;
14846dbb499eSCian Wilson   PetscFunctionReturn(0);
14856dbb499eSCian Wilson }
14866dbb499eSCian Wilson 
14876dbb499eSCian Wilson #undef __FUNCT__
1488704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
14891e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
1490704ba839SBarry Smith {
1491704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1492704ba839SBarry Smith   PetscErrorCode    ierr;
1493704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
1494704ba839SBarry Smith   char              prefix[128];
1495704ba839SBarry Smith 
1496704ba839SBarry Smith   PetscFunctionBegin;
14976c924f48SJed Brown   if (jac->splitdefined) {
14986c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
14996c924f48SJed Brown     PetscFunctionReturn(0);
15006c924f48SJed Brown   }
1501b00a9115SJed Brown   ierr = PetscNew(&ilink);CHKERRQ(ierr);
1502a04f6461SBarry Smith   if (splitname) {
1503db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1504a04f6461SBarry Smith   } else {
1505785e854fSJed Brown     ierr = PetscMalloc1(8,&ilink->splitname);CHKERRQ(ierr);
1506b5787286SJed Brown     ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr);
1507a04f6461SBarry Smith   }
15089df09d43SBarry Smith   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 */
15091ab39975SBarry Smith   ierr          = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1510b5787286SJed Brown   ierr          = ISDestroy(&ilink->is);CHKERRQ(ierr);
1511b5787286SJed Brown   ilink->is     = is;
1512b5787286SJed Brown   ierr          = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1513b5787286SJed Brown   ierr          = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
1514b5787286SJed Brown   ilink->is_col = is;
15150298fd71SBarry Smith   ilink->next   = NULL;
1516ce94432eSBarry Smith   ierr          = KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);CHKERRQ(ierr);
1517422a814eSBarry Smith   ierr          = KSPSetErrorIfNotConverged(ilink->ksp,pc->erroriffailure);CHKERRQ(ierr);
151820252d06SBarry Smith   ierr          = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1519704ba839SBarry Smith   ierr          = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
15209005cf84SBarry Smith   ierr          = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
1521704ba839SBarry Smith 
15228caf3d72SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);CHKERRQ(ierr);
1523704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
1524704ba839SBarry Smith 
1525704ba839SBarry Smith   if (!next) {
1526704ba839SBarry Smith     jac->head       = ilink;
15270298fd71SBarry Smith     ilink->previous = NULL;
1528704ba839SBarry Smith   } else {
1529704ba839SBarry Smith     while (next->next) {
1530704ba839SBarry Smith       next = next->next;
1531704ba839SBarry Smith     }
1532704ba839SBarry Smith     next->next      = ilink;
1533704ba839SBarry Smith     ilink->previous = next;
1534704ba839SBarry Smith   }
1535704ba839SBarry Smith   jac->nsplits++;
1536704ba839SBarry Smith   PetscFunctionReturn(0);
1537704ba839SBarry Smith }
1538704ba839SBarry Smith 
15390971522cSBarry Smith #undef __FUNCT__
15400971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
15410971522cSBarry Smith /*@
15420971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
15430971522cSBarry Smith 
1544ad4df100SBarry Smith     Logically Collective on PC
15450971522cSBarry Smith 
15460971522cSBarry Smith     Input Parameters:
15470971522cSBarry Smith +   pc  - the preconditioner context
15480298fd71SBarry Smith .   splitname - name of this split, if NULL the number of the split is used
15490971522cSBarry Smith .   n - the number of fields in this split
1550db4c96c1SJed Brown -   fields - the fields in this split
15510971522cSBarry Smith 
15520971522cSBarry Smith     Level: intermediate
15530971522cSBarry Smith 
1554d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
1555d32f9abdSBarry Smith 
15567287d2a3SDmitry Karpeev      The PCFieldSplitSetFields() is for defining fields as strided blocks. For example, if the block
1557d32f9abdSBarry Smith      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
1558d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1559d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
1560d32f9abdSBarry Smith 
1561db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
1562db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1563db4c96c1SJed Brown 
15645d4c12cdSJungho Lee      Developer Note: This routine does not actually create the IS representing the split, that is delayed
15655d4c12cdSJungho Lee      until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be
15665d4c12cdSJungho Lee      available when this routine is called.
15675d4c12cdSJungho Lee 
1568d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
15690971522cSBarry Smith 
15700971522cSBarry Smith @*/
15715d4c12cdSJungho Lee PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
15720971522cSBarry Smith {
15734ac538c5SBarry Smith   PetscErrorCode ierr;
15740971522cSBarry Smith 
15750971522cSBarry Smith   PetscFunctionBegin;
15760700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1577db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1578ce94432eSBarry Smith   if (n < 1) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1579db4c96c1SJed Brown   PetscValidIntPointer(fields,3);
15800246f55cSBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt*,const PetscInt*),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr);
15810971522cSBarry Smith   PetscFunctionReturn(0);
15820971522cSBarry Smith }
15830971522cSBarry Smith 
15840971522cSBarry Smith #undef __FUNCT__
1585c84da90fSDmitry Karpeev #define __FUNCT__ "PCFieldSplitSetDiagUseAmat"
1586c84da90fSDmitry Karpeev /*@
1587c84da90fSDmitry Karpeev     PCFieldSplitSetDiagUseAmat - set flag indicating whether to extract diagonal blocks from Amat (rather than Pmat)
1588c84da90fSDmitry Karpeev 
1589c84da90fSDmitry Karpeev     Logically Collective on PC
1590c84da90fSDmitry Karpeev 
1591c84da90fSDmitry Karpeev     Input Parameters:
1592c84da90fSDmitry Karpeev +   pc  - the preconditioner object
1593c84da90fSDmitry Karpeev -   flg - boolean flag indicating whether or not to use Amat to extract the diagonal blocks from
1594c84da90fSDmitry Karpeev 
15959506b623SDmitry Karpeev     Options Database:
15969506b623SDmitry Karpeev .     -pc_fieldsplit_diag_use_amat
1597c84da90fSDmitry Karpeev 
1598c84da90fSDmitry Karpeev     Level: intermediate
1599c84da90fSDmitry Karpeev 
1600c84da90fSDmitry Karpeev .seealso: PCFieldSplitGetDiagUseAmat(), PCFieldSplitSetOffDiagUseAmat(), PCFIELDSPLIT
1601c84da90fSDmitry Karpeev 
1602c84da90fSDmitry Karpeev @*/
1603c84da90fSDmitry Karpeev PetscErrorCode  PCFieldSplitSetDiagUseAmat(PC pc,PetscBool flg)
1604c84da90fSDmitry Karpeev {
1605c84da90fSDmitry Karpeev   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1606c84da90fSDmitry Karpeev   PetscBool      isfs;
1607c84da90fSDmitry Karpeev   PetscErrorCode ierr;
1608c84da90fSDmitry Karpeev 
1609c84da90fSDmitry Karpeev   PetscFunctionBegin;
1610c84da90fSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1611c84da90fSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
1612c84da90fSDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1613c84da90fSDmitry Karpeev   jac->diag_use_amat = flg;
1614c84da90fSDmitry Karpeev   PetscFunctionReturn(0);
1615c84da90fSDmitry Karpeev }
1616c84da90fSDmitry Karpeev 
1617c84da90fSDmitry Karpeev #undef __FUNCT__
1618c84da90fSDmitry Karpeev #define __FUNCT__ "PCFieldSplitGetDiagUseAmat"
1619c84da90fSDmitry Karpeev /*@
1620c84da90fSDmitry Karpeev     PCFieldSplitGetDiagUseAmat - get the flag indicating whether to extract diagonal blocks from Amat (rather than Pmat)
1621c84da90fSDmitry Karpeev 
1622c84da90fSDmitry Karpeev     Logically Collective on PC
1623c84da90fSDmitry Karpeev 
1624c84da90fSDmitry Karpeev     Input Parameters:
1625c84da90fSDmitry Karpeev .   pc  - the preconditioner object
1626c84da90fSDmitry Karpeev 
1627c84da90fSDmitry Karpeev     Output Parameters:
1628c84da90fSDmitry Karpeev .   flg - boolean flag indicating whether or not to use Amat to extract the diagonal blocks from
1629c84da90fSDmitry Karpeev 
1630c84da90fSDmitry Karpeev 
1631c84da90fSDmitry Karpeev     Level: intermediate
1632c84da90fSDmitry Karpeev 
1633c84da90fSDmitry Karpeev .seealso: PCFieldSplitSetDiagUseAmat(), PCFieldSplitGetOffDiagUseAmat(), PCFIELDSPLIT
1634c84da90fSDmitry Karpeev 
1635c84da90fSDmitry Karpeev @*/
1636c84da90fSDmitry Karpeev PetscErrorCode  PCFieldSplitGetDiagUseAmat(PC pc,PetscBool *flg)
1637c84da90fSDmitry Karpeev {
1638c84da90fSDmitry Karpeev   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1639c84da90fSDmitry Karpeev   PetscBool      isfs;
1640c84da90fSDmitry Karpeev   PetscErrorCode ierr;
1641c84da90fSDmitry Karpeev 
1642c84da90fSDmitry Karpeev   PetscFunctionBegin;
1643c84da90fSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1644c84da90fSDmitry Karpeev   PetscValidPointer(flg,2);
1645c84da90fSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
1646c84da90fSDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1647c84da90fSDmitry Karpeev   *flg = jac->diag_use_amat;
1648c84da90fSDmitry Karpeev   PetscFunctionReturn(0);
1649c84da90fSDmitry Karpeev }
1650c84da90fSDmitry Karpeev 
1651c84da90fSDmitry Karpeev #undef __FUNCT__
1652c84da90fSDmitry Karpeev #define __FUNCT__ "PCFieldSplitSetOffDiagUseAmat"
1653c84da90fSDmitry Karpeev /*@
1654c84da90fSDmitry Karpeev     PCFieldSplitSetOffDiagUseAmat - set flag indicating whether to extract off-diagonal blocks from Amat (rather than Pmat)
1655c84da90fSDmitry Karpeev 
1656c84da90fSDmitry Karpeev     Logically Collective on PC
1657c84da90fSDmitry Karpeev 
1658c84da90fSDmitry Karpeev     Input Parameters:
1659c84da90fSDmitry Karpeev +   pc  - the preconditioner object
1660c84da90fSDmitry Karpeev -   flg - boolean flag indicating whether or not to use Amat to extract the off-diagonal blocks from
1661c84da90fSDmitry Karpeev 
16629506b623SDmitry Karpeev     Options Database:
16639506b623SDmitry Karpeev .     -pc_fieldsplit_off_diag_use_amat
1664c84da90fSDmitry Karpeev 
1665c84da90fSDmitry Karpeev     Level: intermediate
1666c84da90fSDmitry Karpeev 
1667c84da90fSDmitry Karpeev .seealso: PCFieldSplitGetOffDiagUseAmat(), PCFieldSplitSetDiagUseAmat(), PCFIELDSPLIT
1668c84da90fSDmitry Karpeev 
1669c84da90fSDmitry Karpeev @*/
1670c84da90fSDmitry Karpeev PetscErrorCode  PCFieldSplitSetOffDiagUseAmat(PC pc,PetscBool flg)
1671c84da90fSDmitry Karpeev {
1672c84da90fSDmitry Karpeev   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1673c84da90fSDmitry Karpeev   PetscBool      isfs;
1674c84da90fSDmitry Karpeev   PetscErrorCode ierr;
1675c84da90fSDmitry Karpeev 
1676c84da90fSDmitry Karpeev   PetscFunctionBegin;
1677c84da90fSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1678c84da90fSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
1679c84da90fSDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1680c84da90fSDmitry Karpeev   jac->offdiag_use_amat = flg;
1681c84da90fSDmitry Karpeev   PetscFunctionReturn(0);
1682c84da90fSDmitry Karpeev }
1683c84da90fSDmitry Karpeev 
1684c84da90fSDmitry Karpeev #undef __FUNCT__
1685c84da90fSDmitry Karpeev #define __FUNCT__ "PCFieldSplitGetOffDiagUseAmat"
1686c84da90fSDmitry Karpeev /*@
1687c84da90fSDmitry Karpeev     PCFieldSplitGetOffDiagUseAmat - get the flag indicating whether to extract off-diagonal blocks from Amat (rather than Pmat)
1688c84da90fSDmitry Karpeev 
1689c84da90fSDmitry Karpeev     Logically Collective on PC
1690c84da90fSDmitry Karpeev 
1691c84da90fSDmitry Karpeev     Input Parameters:
1692c84da90fSDmitry Karpeev .   pc  - the preconditioner object
1693c84da90fSDmitry Karpeev 
1694c84da90fSDmitry Karpeev     Output Parameters:
1695c84da90fSDmitry Karpeev .   flg - boolean flag indicating whether or not to use Amat to extract the off-diagonal blocks from
1696c84da90fSDmitry Karpeev 
1697c84da90fSDmitry Karpeev 
1698c84da90fSDmitry Karpeev     Level: intermediate
1699c84da90fSDmitry Karpeev 
1700c84da90fSDmitry Karpeev .seealso: PCFieldSplitSetOffDiagUseAmat(), PCFieldSplitGetDiagUseAmat(), PCFIELDSPLIT
1701c84da90fSDmitry Karpeev 
1702c84da90fSDmitry Karpeev @*/
1703c84da90fSDmitry Karpeev PetscErrorCode  PCFieldSplitGetOffDiagUseAmat(PC pc,PetscBool *flg)
1704c84da90fSDmitry Karpeev {
1705c84da90fSDmitry Karpeev   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1706c84da90fSDmitry Karpeev   PetscBool      isfs;
1707c84da90fSDmitry Karpeev   PetscErrorCode ierr;
1708c84da90fSDmitry Karpeev 
1709c84da90fSDmitry Karpeev   PetscFunctionBegin;
1710c84da90fSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1711c84da90fSDmitry Karpeev   PetscValidPointer(flg,2);
1712c84da90fSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
1713c84da90fSDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1714c84da90fSDmitry Karpeev   *flg = jac->offdiag_use_amat;
1715c84da90fSDmitry Karpeev   PetscFunctionReturn(0);
1716c84da90fSDmitry Karpeev }
1717c84da90fSDmitry Karpeev 
1718c84da90fSDmitry Karpeev 
1719c84da90fSDmitry Karpeev 
1720c84da90fSDmitry Karpeev #undef __FUNCT__
1721704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
1722218d4943SBarry Smith /*@C
1723704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
1724704ba839SBarry Smith 
1725ad4df100SBarry Smith     Logically Collective on PC
1726704ba839SBarry Smith 
1727704ba839SBarry Smith     Input Parameters:
1728704ba839SBarry Smith +   pc  - the preconditioner context
17290298fd71SBarry Smith .   splitname - name of this split, if NULL the number of the split is used
1730db4c96c1SJed Brown -   is - the index set that defines the vector elements in this field
1731704ba839SBarry Smith 
1732d32f9abdSBarry Smith 
1733a6ffb8dbSJed Brown     Notes:
1734a6ffb8dbSJed Brown     Use PCFieldSplitSetFields(), for fields defined by strided types.
1735a6ffb8dbSJed Brown 
1736db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
1737db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1738d32f9abdSBarry Smith 
1739704ba839SBarry Smith     Level: intermediate
1740704ba839SBarry Smith 
1741704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
1742704ba839SBarry Smith 
1743704ba839SBarry Smith @*/
17447087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1745704ba839SBarry Smith {
17464ac538c5SBarry Smith   PetscErrorCode ierr;
1747704ba839SBarry Smith 
1748704ba839SBarry Smith   PetscFunctionBegin;
17490700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
17507b62db95SJungho Lee   if (splitname) PetscValidCharPointer(splitname,2);
1751db4c96c1SJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,3);
1752ccc738f7SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr);
1753704ba839SBarry Smith   PetscFunctionReturn(0);
1754704ba839SBarry Smith }
1755704ba839SBarry Smith 
1756704ba839SBarry Smith #undef __FUNCT__
175757a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS"
175857a9adfeSMatthew G Knepley /*@
175957a9adfeSMatthew G Knepley     PCFieldSplitGetIS - Retrieves the elements for a field as an IS
176057a9adfeSMatthew G Knepley 
176157a9adfeSMatthew G Knepley     Logically Collective on PC
176257a9adfeSMatthew G Knepley 
176357a9adfeSMatthew G Knepley     Input Parameters:
176457a9adfeSMatthew G Knepley +   pc  - the preconditioner context
176557a9adfeSMatthew G Knepley -   splitname - name of this split
176657a9adfeSMatthew G Knepley 
176757a9adfeSMatthew G Knepley     Output Parameter:
17680298fd71SBarry Smith -   is - the index set that defines the vector elements in this field, or NULL if the field is not found
176957a9adfeSMatthew G Knepley 
177057a9adfeSMatthew G Knepley     Level: intermediate
177157a9adfeSMatthew G Knepley 
177257a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS()
177357a9adfeSMatthew G Knepley 
177457a9adfeSMatthew G Knepley @*/
177557a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
177657a9adfeSMatthew G Knepley {
177757a9adfeSMatthew G Knepley   PetscErrorCode ierr;
177857a9adfeSMatthew G Knepley 
177957a9adfeSMatthew G Knepley   PetscFunctionBegin;
178057a9adfeSMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
178157a9adfeSMatthew G Knepley   PetscValidCharPointer(splitname,2);
178257a9adfeSMatthew G Knepley   PetscValidPointer(is,3);
178357a9adfeSMatthew G Knepley   {
178457a9adfeSMatthew G Knepley     PC_FieldSplit     *jac  = (PC_FieldSplit*) pc->data;
178557a9adfeSMatthew G Knepley     PC_FieldSplitLink ilink = jac->head;
178657a9adfeSMatthew G Knepley     PetscBool         found;
178757a9adfeSMatthew G Knepley 
17880298fd71SBarry Smith     *is = NULL;
178957a9adfeSMatthew G Knepley     while (ilink) {
179057a9adfeSMatthew G Knepley       ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr);
179157a9adfeSMatthew G Knepley       if (found) {
179257a9adfeSMatthew G Knepley         *is = ilink->is;
179357a9adfeSMatthew G Knepley         break;
179457a9adfeSMatthew G Knepley       }
179557a9adfeSMatthew G Knepley       ilink = ilink->next;
179657a9adfeSMatthew G Knepley     }
179757a9adfeSMatthew G Knepley   }
179857a9adfeSMatthew G Knepley   PetscFunctionReturn(0);
179957a9adfeSMatthew G Knepley }
180057a9adfeSMatthew G Knepley 
180157a9adfeSMatthew G Knepley #undef __FUNCT__
180251f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
180351f519a2SBarry Smith /*@
180451f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
180551f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
180651f519a2SBarry Smith 
1807ad4df100SBarry Smith     Logically Collective on PC
180851f519a2SBarry Smith 
180951f519a2SBarry Smith     Input Parameters:
181051f519a2SBarry Smith +   pc  - the preconditioner context
181151f519a2SBarry Smith -   bs - the block size
181251f519a2SBarry Smith 
181351f519a2SBarry Smith     Level: intermediate
181451f519a2SBarry Smith 
181551f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
181651f519a2SBarry Smith 
181751f519a2SBarry Smith @*/
18187087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
181951f519a2SBarry Smith {
18204ac538c5SBarry Smith   PetscErrorCode ierr;
182151f519a2SBarry Smith 
182251f519a2SBarry Smith   PetscFunctionBegin;
18230700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1824c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,bs,2);
18254ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr);
182651f519a2SBarry Smith   PetscFunctionReturn(0);
182751f519a2SBarry Smith }
182851f519a2SBarry Smith 
182951f519a2SBarry Smith #undef __FUNCT__
183069a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
18310971522cSBarry Smith /*@C
183269a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
18330971522cSBarry Smith 
183469a612a9SBarry Smith    Collective on KSP
18350971522cSBarry Smith 
18360971522cSBarry Smith    Input Parameter:
18370971522cSBarry Smith .  pc - the preconditioner context
18380971522cSBarry Smith 
18390971522cSBarry Smith    Output Parameters:
184013e0d083SBarry Smith +  n - the number of splits
184169a612a9SBarry Smith -  pc - the array of KSP contexts
18420971522cSBarry Smith 
18430971522cSBarry Smith    Note:
1844d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
1845d32f9abdSBarry Smith    (not the KSP just the array that contains them).
18460971522cSBarry Smith 
184769a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
18480971522cSBarry Smith 
1849196cc216SBarry Smith    Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs.
18500298fd71SBarry Smith       You can call PCFieldSplitGetSubKSP(pc,n,NULL_OBJECT,ierr) to determine how large the
1851196cc216SBarry Smith       KSP array must be.
1852196cc216SBarry Smith 
1853196cc216SBarry Smith 
18540971522cSBarry Smith    Level: advanced
18550971522cSBarry Smith 
18560971522cSBarry Smith .seealso: PCFIELDSPLIT
18570971522cSBarry Smith @*/
18587087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
18590971522cSBarry Smith {
18604ac538c5SBarry Smith   PetscErrorCode ierr;
18610971522cSBarry Smith 
18620971522cSBarry Smith   PetscFunctionBegin;
18630700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
186413e0d083SBarry Smith   if (n) PetscValidIntPointer(n,2);
18654ac538c5SBarry Smith   ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
18660971522cSBarry Smith   PetscFunctionReturn(0);
18670971522cSBarry Smith }
18680971522cSBarry Smith 
1869e69d4d44SBarry Smith #undef __FUNCT__
187029f8a81cSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurPre"
1871e69d4d44SBarry Smith /*@
187229f8a81cSJed Brown     PCFieldSplitSetSchurPre -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
1873a04f6461SBarry Smith       A11 matrix. Otherwise no preconditioner is used.
1874e69d4d44SBarry Smith 
1875e69d4d44SBarry Smith     Collective on PC
1876e69d4d44SBarry Smith 
1877e69d4d44SBarry Smith     Input Parameters:
1878e69d4d44SBarry Smith +   pc      - the preconditioner context
18796fdd48a9SBarry Smith .   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
18806fdd48a9SBarry Smith               PC_FIELDSPLIT_SCHUR_PRE_SELFP, and PC_FIELDSPLIT_SCHUR_PRE_FULL
18810298fd71SBarry Smith -   userpre - matrix to use for preconditioning, or NULL
1882084e4875SJed Brown 
1883e69d4d44SBarry Smith     Options Database:
18846fdd48a9SBarry Smith .     -pc_fieldsplit_schur_precondition <self,selfp,user,a11,full> - default is a11. See notes for meaning of various arguments
1885e69d4d44SBarry Smith 
1886fd1303e9SJungho Lee     Notes:
1887fd1303e9SJungho Lee $    If ptype is
1888a6a584a2SBarry Smith $        a11 then the preconditioner for the Schur complement is generated from the block diagonal part of the preconditioner
1889585af1b9SMatthew G. Knepley $             matrix associated with the Schur complement (i.e. A11), not he Schur complement matrix
1890a6a584a2SBarry Smith $        self the preconditioner for the Schur complement is generated from the symbolic representation of the Schur complement matrix:
1891a6a584a2SBarry Smith $             The only preconditioner that currently works with this symbolic respresentation matrix object is the PCLSC
1892fd1303e9SJungho Lee $             preconditioner
18936fdd48a9SBarry Smith $        user then the preconditioner for the Schur complement is generated from the user provided matrix (pre argument
18946fdd48a9SBarry Smith $             to this function).
1895a6a584a2SBarry Smith $        selfp then the preconditioning for the Schur complement is generated from an explicitly-assembled approximation Sp = A11 - A10 inv(diag(A00)) A01
18961d71051fSDmitry Karpeev $             This is only a good preconditioner when diag(A00) is a good preconditioner for A00. Optionally, A00 can be
18976fdd48a9SBarry Smith $             lumped before extracting the diagonal using the additional option -fieldsplit_1_mat_schur_complement_ainv_type lump
18986fdd48a9SBarry Smith $        full then the preconditioner for the Schur complement is generated from the exact Schur complement matrix representation computed internally by PFIELDSPLIT (this is expensive)
18996fdd48a9SBarry Smith $             useful mostly as a test that the Schur complement approach can work for your problem
1900fd1303e9SJungho Lee 
1901e87fab1bSBarry Smith      When solving a saddle point problem, where the A11 block is identically zero, using a11 as the ptype only makes sense
1902fd1303e9SJungho Lee     with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
1903a7476a74SDmitry Karpeev     -fieldsplit_1_pc_type lsc which uses the least squares commutator to compute a preconditioner for the Schur complement.
1904fd1303e9SJungho Lee 
1905e69d4d44SBarry Smith     Level: intermediate
1906e69d4d44SBarry Smith 
19071d71051fSDmitry Karpeev .seealso: PCFieldSplitGetSchurPre(), PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType,
19081d71051fSDmitry Karpeev           MatSchurComplementSetAinvType(), PCLSC
1909e69d4d44SBarry Smith 
1910e69d4d44SBarry Smith @*/
191129f8a81cSJed Brown PetscErrorCode PCFieldSplitSetSchurPre(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1912e69d4d44SBarry Smith {
19134ac538c5SBarry Smith   PetscErrorCode ierr;
1914e69d4d44SBarry Smith 
1915e69d4d44SBarry Smith   PetscFunctionBegin;
19160700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
191729f8a81cSJed Brown   ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurPre_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr);
1918e69d4d44SBarry Smith   PetscFunctionReturn(0);
1919e69d4d44SBarry Smith }
192029f8a81cSJed Brown PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) {return PCFieldSplitSetSchurPre(pc,ptype,pre);} /* Deprecated name */
1921e69d4d44SBarry Smith 
1922e69d4d44SBarry Smith #undef __FUNCT__
192337a82bf0SJed Brown #define __FUNCT__ "PCFieldSplitGetSchurPre"
192437a82bf0SJed Brown /*@
192537a82bf0SJed Brown     PCFieldSplitGetSchurPre - For Schur complement fieldsplit, determine how the Schur complement will be
192629f8a81cSJed Brown     preconditioned.  See PCFieldSplitSetSchurPre() for details.
192737a82bf0SJed Brown 
192837a82bf0SJed Brown     Logically Collective on PC
192937a82bf0SJed Brown 
193037a82bf0SJed Brown     Input Parameters:
193137a82bf0SJed Brown .   pc      - the preconditioner context
193237a82bf0SJed Brown 
193337a82bf0SJed Brown     Output Parameters:
193437a82bf0SJed Brown +   ptype   - which matrix to use for preconditioning the Schur complement: PC_FIELDSPLIT_SCHUR_PRE_A11, PC_FIELDSPLIT_SCHUR_PRE_SELF, PC_FIELDSPLIT_PRE_USER
193537a82bf0SJed Brown -   userpre - matrix to use for preconditioning (with PC_FIELDSPLIT_PRE_USER), or NULL
193637a82bf0SJed Brown 
193737a82bf0SJed Brown     Level: intermediate
193837a82bf0SJed Brown 
193929f8a81cSJed Brown .seealso: PCFieldSplitSetSchurPre(), PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC
194037a82bf0SJed Brown 
194137a82bf0SJed Brown @*/
194237a82bf0SJed Brown PetscErrorCode PCFieldSplitGetSchurPre(PC pc,PCFieldSplitSchurPreType *ptype,Mat *pre)
194337a82bf0SJed Brown {
194437a82bf0SJed Brown   PetscErrorCode ierr;
194537a82bf0SJed Brown 
194637a82bf0SJed Brown   PetscFunctionBegin;
194737a82bf0SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
194837a82bf0SJed Brown   ierr = PetscUseMethod(pc,"PCFieldSplitGetSchurPre_C",(PC,PCFieldSplitSchurPreType*,Mat*),(pc,ptype,pre));CHKERRQ(ierr);
1949e69d4d44SBarry Smith   PetscFunctionReturn(0);
1950e69d4d44SBarry Smith }
1951e69d4d44SBarry Smith 
1952e69d4d44SBarry Smith #undef __FUNCT__
1953470b340bSDmitry Karpeev #define __FUNCT__ "PCFieldSplitSchurGetS"
195445e7fc46SDmitry Karpeev /*@
1955470b340bSDmitry Karpeev     PCFieldSplitSchurGetS -  extract the MatSchurComplement object used by this PC in case it needs to be configured separately
195645e7fc46SDmitry Karpeev 
195745e7fc46SDmitry Karpeev     Not collective
195845e7fc46SDmitry Karpeev 
195945e7fc46SDmitry Karpeev     Input Parameter:
196045e7fc46SDmitry Karpeev .   pc      - the preconditioner context
196145e7fc46SDmitry Karpeev 
196245e7fc46SDmitry Karpeev     Output Parameter:
1963470b340bSDmitry Karpeev .   S       - the Schur complement matrix
196445e7fc46SDmitry Karpeev 
1965470b340bSDmitry Karpeev     Notes:
1966470b340bSDmitry Karpeev     This matrix should not be destroyed using MatDestroy(); rather, use PCFieldSplitSchurRestoreS().
196745e7fc46SDmitry Karpeev 
196845e7fc46SDmitry Karpeev     Level: advanced
196945e7fc46SDmitry Karpeev 
197029f8a81cSJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSchurPreType, PCFieldSplitSetSchurPre(), MatSchurComplement, PCFieldSplitSchurRestoreS()
197145e7fc46SDmitry Karpeev 
197245e7fc46SDmitry Karpeev @*/
1973470b340bSDmitry Karpeev PetscErrorCode  PCFieldSplitSchurGetS(PC pc,Mat *S)
197445e7fc46SDmitry Karpeev {
197545e7fc46SDmitry Karpeev   PetscErrorCode ierr;
197645e7fc46SDmitry Karpeev   const char*    t;
197745e7fc46SDmitry Karpeev   PetscBool      isfs;
197845e7fc46SDmitry Karpeev   PC_FieldSplit  *jac;
197945e7fc46SDmitry Karpeev 
198045e7fc46SDmitry Karpeev   PetscFunctionBegin;
198145e7fc46SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
198245e7fc46SDmitry Karpeev   ierr = PetscObjectGetType((PetscObject)pc,&t);CHKERRQ(ierr);
198345e7fc46SDmitry Karpeev   ierr = PetscStrcmp(t,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
198445e7fc46SDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PC of type PCFIELDSPLIT, got %s instead",t);
198545e7fc46SDmitry Karpeev   jac = (PC_FieldSplit*)pc->data;
1986470b340bSDmitry Karpeev   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PCFIELDSPLIT of type SCHUR, got %D instead",jac->type);
1987470b340bSDmitry Karpeev   if (S) *S = jac->schur;
198845e7fc46SDmitry Karpeev   PetscFunctionReturn(0);
198945e7fc46SDmitry Karpeev }
199045e7fc46SDmitry Karpeev 
199145e7fc46SDmitry Karpeev #undef __FUNCT__
1992470b340bSDmitry Karpeev #define __FUNCT__ "PCFieldSplitSchurRestoreS"
1993470b340bSDmitry Karpeev /*@
1994470b340bSDmitry Karpeev     PCFieldSplitSchurRestoreS -  restores the MatSchurComplement object used by this PC
1995470b340bSDmitry Karpeev 
1996470b340bSDmitry Karpeev     Not collective
1997470b340bSDmitry Karpeev 
1998470b340bSDmitry Karpeev     Input Parameters:
1999470b340bSDmitry Karpeev +   pc      - the preconditioner context
2000470b340bSDmitry Karpeev .   S       - the Schur complement matrix
2001470b340bSDmitry Karpeev 
2002470b340bSDmitry Karpeev     Level: advanced
2003470b340bSDmitry Karpeev 
200429f8a81cSJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSchurPreType, PCFieldSplitSetSchurPre(), MatSchurComplement, PCFieldSplitSchurGetS()
2005470b340bSDmitry Karpeev 
2006470b340bSDmitry Karpeev @*/
2007bca69d2bSDmitry Karpeev PetscErrorCode  PCFieldSplitSchurRestoreS(PC pc,Mat *S)
2008470b340bSDmitry Karpeev {
2009470b340bSDmitry Karpeev   PetscErrorCode ierr;
2010470b340bSDmitry Karpeev   const char*    t;
2011470b340bSDmitry Karpeev   PetscBool      isfs;
2012470b340bSDmitry Karpeev   PC_FieldSplit  *jac;
2013470b340bSDmitry Karpeev 
2014470b340bSDmitry Karpeev   PetscFunctionBegin;
2015470b340bSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2016470b340bSDmitry Karpeev   ierr = PetscObjectGetType((PetscObject)pc,&t);CHKERRQ(ierr);
2017470b340bSDmitry Karpeev   ierr = PetscStrcmp(t,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
2018470b340bSDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PC of type PCFIELDSPLIT, got %s instead",t);
2019470b340bSDmitry Karpeev   jac = (PC_FieldSplit*)pc->data;
2020470b340bSDmitry Karpeev   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PCFIELDSPLIT of type SCHUR, got %D instead",jac->type);
2021bca69d2bSDmitry Karpeev   if (!S || *S != jac->schur) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"MatSchurComplement restored is not the same as gotten");
2022470b340bSDmitry Karpeev   PetscFunctionReturn(0);
2023470b340bSDmitry Karpeev }
2024470b340bSDmitry Karpeev 
2025470b340bSDmitry Karpeev 
2026470b340bSDmitry Karpeev #undef __FUNCT__
202729f8a81cSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurPre_FieldSplit"
202829f8a81cSJed Brown static PetscErrorCode  PCFieldSplitSetSchurPre_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
2029e69d4d44SBarry Smith {
2030e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
2031084e4875SJed Brown   PetscErrorCode ierr;
2032e69d4d44SBarry Smith 
2033e69d4d44SBarry Smith   PetscFunctionBegin;
2034084e4875SJed Brown   jac->schurpre = ptype;
2035a7476a74SDmitry Karpeev   if (ptype == PC_FIELDSPLIT_SCHUR_PRE_USER && pre) {
20366bf464f9SBarry Smith     ierr            = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
2037084e4875SJed Brown     jac->schur_user = pre;
2038084e4875SJed Brown     ierr            = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
2039084e4875SJed Brown   }
2040e69d4d44SBarry Smith   PetscFunctionReturn(0);
2041e69d4d44SBarry Smith }
2042e69d4d44SBarry Smith 
204330ad9308SMatthew Knepley #undef __FUNCT__
204437a82bf0SJed Brown #define __FUNCT__ "PCFieldSplitGetSchurPre_FieldSplit"
204537a82bf0SJed Brown static PetscErrorCode  PCFieldSplitGetSchurPre_FieldSplit(PC pc,PCFieldSplitSchurPreType *ptype,Mat *pre)
204637a82bf0SJed Brown {
204737a82bf0SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
204837a82bf0SJed Brown 
204937a82bf0SJed Brown   PetscFunctionBegin;
205037a82bf0SJed Brown   *ptype = jac->schurpre;
205137a82bf0SJed Brown   *pre   = jac->schur_user;
205237a82bf0SJed Brown   PetscFunctionReturn(0);
205337a82bf0SJed Brown }
205437a82bf0SJed Brown 
205537a82bf0SJed Brown #undef __FUNCT__
2056c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType"
2057ab1df9f5SJed Brown /*@
2058c9c6ffaaSJed Brown     PCFieldSplitSetSchurFactType -  sets which blocks of the approximate block factorization to retain
2059ab1df9f5SJed Brown 
2060ab1df9f5SJed Brown     Collective on PC
2061ab1df9f5SJed Brown 
2062ab1df9f5SJed Brown     Input Parameters:
2063ab1df9f5SJed Brown +   pc  - the preconditioner context
2064c9c6ffaaSJed Brown -   ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default
2065ab1df9f5SJed Brown 
2066ab1df9f5SJed Brown     Options Database:
2067c9c6ffaaSJed Brown .     -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full
2068ab1df9f5SJed Brown 
2069ab1df9f5SJed Brown 
2070ab1df9f5SJed Brown     Level: intermediate
2071ab1df9f5SJed Brown 
2072ab1df9f5SJed Brown     Notes:
2073ab1df9f5SJed Brown     The FULL factorization is
2074ab1df9f5SJed Brown 
2075ab1df9f5SJed Brown $   (A   B)  = (1       0) (A   0) (1  Ainv*B)
2076ab1df9f5SJed Brown $   (C   D)    (C*Ainv  1) (0   S) (0     1  )
2077ab1df9f5SJed Brown 
20786be4592eSBarry Smith     where S = D - 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,
20796be4592eSBarry Smith     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).
2080ab1df9f5SJed Brown 
20816be4592eSBarry Smith     If applied exactly, FULL factorization is a direct solver. The preconditioned operator with LOWER or UPPER has all eigenvalues equal to 1 and minimal polynomial
20826be4592eSBarry Smith     of degree 2, so KSPGMRES converges in 2 iterations. If the iteration count is very low, consider using KSPFGMRES or KSPGCR which can use one less preconditioner
20836be4592eSBarry Smith     application in this case. Note that the preconditioned operator may be highly non-normal, so such fast convergence may not be observed in practice. With DIAG,
20846be4592eSBarry Smith     the preconditioned operator has three distinct nonzero eigenvalues and minimal polynomial of degree at most 4, so KSPGMRES converges in at most 4 iterations.
2085ab1df9f5SJed Brown 
20866be4592eSBarry Smith     For symmetric problems in which A is positive definite and S is negative definite, DIAG can be used with KSPMINRES. Note that a flexible method like KSPFGMRES
20876be4592eSBarry Smith     or KSPGCR must be used if the fieldsplit preconditioner is nonlinear (e.g. a few iterations of a Krylov method is used inside a split).
2088ab1df9f5SJed Brown 
2089ab1df9f5SJed Brown     References:
209096a0c994SBarry Smith +   1. - Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000).
209196a0c994SBarry Smith -   2. - Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001).
2092ab1df9f5SJed Brown 
2093ab1df9f5SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType
2094ab1df9f5SJed Brown @*/
2095c9c6ffaaSJed Brown PetscErrorCode  PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype)
2096ab1df9f5SJed Brown {
2097ab1df9f5SJed Brown   PetscErrorCode ierr;
2098ab1df9f5SJed Brown 
2099ab1df9f5SJed Brown   PetscFunctionBegin;
2100ab1df9f5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2101c9c6ffaaSJed Brown   ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr);
2102ab1df9f5SJed Brown   PetscFunctionReturn(0);
2103ab1df9f5SJed Brown }
2104ab1df9f5SJed Brown 
2105ab1df9f5SJed Brown #undef __FUNCT__
2106c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType_FieldSplit"
21071e6b0712SBarry Smith static PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype)
2108ab1df9f5SJed Brown {
2109ab1df9f5SJed Brown   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2110ab1df9f5SJed Brown 
2111ab1df9f5SJed Brown   PetscFunctionBegin;
2112ab1df9f5SJed Brown   jac->schurfactorization = ftype;
2113ab1df9f5SJed Brown   PetscFunctionReturn(0);
2114ab1df9f5SJed Brown }
2115ab1df9f5SJed Brown 
2116ab1df9f5SJed Brown #undef __FUNCT__
211730ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
211830ad9308SMatthew Knepley /*@C
21198c03b21aSDmitry Karpeev    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
212030ad9308SMatthew Knepley 
212130ad9308SMatthew Knepley    Collective on KSP
212230ad9308SMatthew Knepley 
212330ad9308SMatthew Knepley    Input Parameter:
212430ad9308SMatthew Knepley .  pc - the preconditioner context
212530ad9308SMatthew Knepley 
212630ad9308SMatthew Knepley    Output Parameters:
2127a04f6461SBarry Smith +  A00 - the (0,0) block
2128a04f6461SBarry Smith .  A01 - the (0,1) block
2129a04f6461SBarry Smith .  A10 - the (1,0) block
2130a04f6461SBarry Smith -  A11 - the (1,1) block
213130ad9308SMatthew Knepley 
213230ad9308SMatthew Knepley    Level: advanced
213330ad9308SMatthew Knepley 
213430ad9308SMatthew Knepley .seealso: PCFIELDSPLIT
213530ad9308SMatthew Knepley @*/
2136a04f6461SBarry Smith PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
213730ad9308SMatthew Knepley {
213830ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;
213930ad9308SMatthew Knepley 
214030ad9308SMatthew Knepley   PetscFunctionBegin;
21410700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2142ce94432eSBarry Smith   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
2143a04f6461SBarry Smith   if (A00) *A00 = jac->pmat[0];
2144a04f6461SBarry Smith   if (A01) *A01 = jac->B;
2145a04f6461SBarry Smith   if (A10) *A10 = jac->C;
2146a04f6461SBarry Smith   if (A11) *A11 = jac->pmat[1];
214730ad9308SMatthew Knepley   PetscFunctionReturn(0);
214830ad9308SMatthew Knepley }
214930ad9308SMatthew Knepley 
215079416396SBarry Smith #undef __FUNCT__
215179416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
21521e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
215379416396SBarry Smith {
215479416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
2155e69d4d44SBarry Smith   PetscErrorCode ierr;
215679416396SBarry Smith 
215779416396SBarry Smith   PetscFunctionBegin;
215879416396SBarry Smith   jac->type = type;
21593b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
21603b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
21613b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
21622fa5cd67SKarl Rupp 
2163bdf89e91SBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
216429f8a81cSJed Brown     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",PCFieldSplitSetSchurPre_FieldSplit);CHKERRQ(ierr);
216537a82bf0SJed Brown     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",PCFieldSplitGetSchurPre_FieldSplit);CHKERRQ(ierr);
2166bdf89e91SBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr);
2167e69d4d44SBarry Smith 
21683b224e63SBarry Smith   } else {
21693b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
21703b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
21712fa5cd67SKarl Rupp 
2172bdf89e91SBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
217329f8a81cSJed Brown     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",0);CHKERRQ(ierr);
217437a82bf0SJed Brown     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",0);CHKERRQ(ierr);
2175bdf89e91SBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",0);CHKERRQ(ierr);
21763b224e63SBarry Smith   }
217779416396SBarry Smith   PetscFunctionReturn(0);
217879416396SBarry Smith }
217979416396SBarry Smith 
218051f519a2SBarry Smith #undef __FUNCT__
218151f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
21821e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
218351f519a2SBarry Smith {
218451f519a2SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
218551f519a2SBarry Smith 
218651f519a2SBarry Smith   PetscFunctionBegin;
2187ce94432eSBarry Smith   if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
2188ce94432eSBarry Smith   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);
218951f519a2SBarry Smith   jac->bs = bs;
219051f519a2SBarry Smith   PetscFunctionReturn(0);
219151f519a2SBarry Smith }
219251f519a2SBarry Smith 
219379416396SBarry Smith #undef __FUNCT__
219479416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
2195bc08b0f1SBarry Smith /*@
219679416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
219779416396SBarry Smith 
219879416396SBarry Smith    Collective on PC
219979416396SBarry Smith 
220079416396SBarry Smith    Input Parameter:
220179416396SBarry Smith .  pc - the preconditioner context
220281540f2fSBarry Smith .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
220379416396SBarry Smith 
220479416396SBarry Smith    Options Database Key:
2205a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
220679416396SBarry Smith 
2207b02e2d75SMatthew G Knepley    Level: Intermediate
220879416396SBarry Smith 
220979416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
221079416396SBarry Smith 
221179416396SBarry Smith .seealso: PCCompositeSetType()
221279416396SBarry Smith 
221379416396SBarry Smith @*/
22147087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
221579416396SBarry Smith {
22164ac538c5SBarry Smith   PetscErrorCode ierr;
221779416396SBarry Smith 
221879416396SBarry Smith   PetscFunctionBegin;
22190700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
22204ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr);
222179416396SBarry Smith   PetscFunctionReturn(0);
222279416396SBarry Smith }
222379416396SBarry Smith 
2224b02e2d75SMatthew G Knepley #undef __FUNCT__
2225b02e2d75SMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetType"
2226b02e2d75SMatthew G Knepley /*@
2227b02e2d75SMatthew G Knepley   PCFieldSplitGetType - Gets the type of fieldsplit preconditioner.
2228b02e2d75SMatthew G Knepley 
2229b02e2d75SMatthew G Knepley   Not collective
2230b02e2d75SMatthew G Knepley 
2231b02e2d75SMatthew G Knepley   Input Parameter:
2232b02e2d75SMatthew G Knepley . pc - the preconditioner context
2233b02e2d75SMatthew G Knepley 
2234b02e2d75SMatthew G Knepley   Output Parameter:
2235b02e2d75SMatthew G Knepley . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
2236b02e2d75SMatthew G Knepley 
2237b02e2d75SMatthew G Knepley   Level: Intermediate
2238b02e2d75SMatthew G Knepley 
2239b02e2d75SMatthew G Knepley .keywords: PC, set, type, composite preconditioner, additive, multiplicative
2240b02e2d75SMatthew G Knepley .seealso: PCCompositeSetType()
2241b02e2d75SMatthew G Knepley @*/
2242b02e2d75SMatthew G Knepley PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type)
2243b02e2d75SMatthew G Knepley {
2244b02e2d75SMatthew G Knepley   PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;
2245b02e2d75SMatthew G Knepley 
2246b02e2d75SMatthew G Knepley   PetscFunctionBegin;
2247b02e2d75SMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2248b02e2d75SMatthew G Knepley   PetscValidIntPointer(type,2);
2249b02e2d75SMatthew G Knepley   *type = jac->type;
2250b02e2d75SMatthew G Knepley   PetscFunctionReturn(0);
2251b02e2d75SMatthew G Knepley }
2252b02e2d75SMatthew G Knepley 
22534ab8060aSDmitry Karpeev #undef __FUNCT__
22544ab8060aSDmitry Karpeev #define __FUNCT__ "PCFieldSplitSetDMSplits"
22554ab8060aSDmitry Karpeev /*@
22564ab8060aSDmitry Karpeev    PCFieldSplitSetDMSplits - Flags whether DMCreateFieldDecomposition() should be used to define the splits, whenever possible.
22574ab8060aSDmitry Karpeev 
22584ab8060aSDmitry Karpeev    Logically Collective
22594ab8060aSDmitry Karpeev 
22604ab8060aSDmitry Karpeev    Input Parameters:
22614ab8060aSDmitry Karpeev +  pc   - the preconditioner context
22624ab8060aSDmitry Karpeev -  flg  - boolean indicating whether to use field splits defined by the DM
22634ab8060aSDmitry Karpeev 
22644ab8060aSDmitry Karpeev    Options Database Key:
22654ab8060aSDmitry Karpeev .  -pc_fieldsplit_dm_splits
22664ab8060aSDmitry Karpeev 
22674ab8060aSDmitry Karpeev    Level: Intermediate
22684ab8060aSDmitry Karpeev 
22694ab8060aSDmitry Karpeev .keywords: PC, DM, composite preconditioner, additive, multiplicative
22704ab8060aSDmitry Karpeev 
22714ab8060aSDmitry Karpeev .seealso: PCFieldSplitGetDMSplits()
22724ab8060aSDmitry Karpeev 
22734ab8060aSDmitry Karpeev @*/
22744ab8060aSDmitry Karpeev PetscErrorCode  PCFieldSplitSetDMSplits(PC pc,PetscBool flg)
22754ab8060aSDmitry Karpeev {
22764ab8060aSDmitry Karpeev   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
22774ab8060aSDmitry Karpeev   PetscBool      isfs;
22784ab8060aSDmitry Karpeev   PetscErrorCode ierr;
22794ab8060aSDmitry Karpeev 
22804ab8060aSDmitry Karpeev   PetscFunctionBegin;
22814ab8060aSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
22824ab8060aSDmitry Karpeev   PetscValidLogicalCollectiveBool(pc,flg,2);
22834ab8060aSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
22844ab8060aSDmitry Karpeev   if (isfs) {
22854ab8060aSDmitry Karpeev     jac->dm_splits = flg;
22864ab8060aSDmitry Karpeev   }
22874ab8060aSDmitry Karpeev   PetscFunctionReturn(0);
22884ab8060aSDmitry Karpeev }
22894ab8060aSDmitry Karpeev 
22904ab8060aSDmitry Karpeev 
22914ab8060aSDmitry Karpeev #undef __FUNCT__
22924ab8060aSDmitry Karpeev #define __FUNCT__ "PCFieldSplitGetDMSplits"
22934ab8060aSDmitry Karpeev /*@
22944ab8060aSDmitry Karpeev    PCFieldSplitGetDMSplits - Returns flag indicating whether DMCreateFieldDecomposition() should be used to define the splits, whenever possible.
22954ab8060aSDmitry Karpeev 
22964ab8060aSDmitry Karpeev    Logically Collective
22974ab8060aSDmitry Karpeev 
22984ab8060aSDmitry Karpeev    Input Parameter:
22994ab8060aSDmitry Karpeev .  pc   - the preconditioner context
23004ab8060aSDmitry Karpeev 
23014ab8060aSDmitry Karpeev    Output Parameter:
23024ab8060aSDmitry Karpeev .  flg  - boolean indicating whether to use field splits defined by the DM
23034ab8060aSDmitry Karpeev 
23044ab8060aSDmitry Karpeev    Level: Intermediate
23054ab8060aSDmitry Karpeev 
23064ab8060aSDmitry Karpeev .keywords: PC, DM, composite preconditioner, additive, multiplicative
23074ab8060aSDmitry Karpeev 
23084ab8060aSDmitry Karpeev .seealso: PCFieldSplitSetDMSplits()
23094ab8060aSDmitry Karpeev 
23104ab8060aSDmitry Karpeev @*/
23114ab8060aSDmitry Karpeev PetscErrorCode  PCFieldSplitGetDMSplits(PC pc,PetscBool* flg)
23124ab8060aSDmitry Karpeev {
23134ab8060aSDmitry Karpeev   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
23144ab8060aSDmitry Karpeev   PetscBool      isfs;
23154ab8060aSDmitry Karpeev   PetscErrorCode ierr;
23164ab8060aSDmitry Karpeev 
23174ab8060aSDmitry Karpeev   PetscFunctionBegin;
23184ab8060aSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
23194ab8060aSDmitry Karpeev   PetscValidPointer(flg,2);
23204ab8060aSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
23214ab8060aSDmitry Karpeev   if (isfs) {
23224ab8060aSDmitry Karpeev     if(flg) *flg = jac->dm_splits;
23234ab8060aSDmitry Karpeev   }
23244ab8060aSDmitry Karpeev   PetscFunctionReturn(0);
23254ab8060aSDmitry Karpeev }
23264ab8060aSDmitry Karpeev 
23270971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
23280971522cSBarry Smith /*MC
2329a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
2330a04f6461SBarry Smith                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.
23310971522cSBarry Smith 
2332edf189efSBarry Smith      To set options on the solvers for each block append -fieldsplit_ to all the PC
2333edf189efSBarry Smith         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
23340971522cSBarry Smith 
2335a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
233669a612a9SBarry Smith          and set the options directly on the resulting KSP object
23370971522cSBarry Smith 
23380971522cSBarry Smith    Level: intermediate
23390971522cSBarry Smith 
234079416396SBarry Smith    Options Database Keys:
234181540f2fSBarry Smith +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
234281540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
234381540f2fSBarry Smith                               been supplied explicitly by -pc_fieldsplit_%d_fields
234481540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
23450f188ba9SJed Brown .   -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur> - type of relaxation or factorization splitting
2346a6a584a2SBarry Smith .   -pc_fieldsplit_schur_precondition <self,selfp,user,a11,full> - default is a11; see PCFieldSplitSetSchurPre()
2347435f959eSBarry Smith .   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver
234879416396SBarry Smith 
23495d4c12cdSJungho Lee -    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
23505d4c12cdSJungho Lee      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
23515d4c12cdSJungho Lee 
2352c8a0d604SMatthew G Knepley    Notes:
2353c8a0d604SMatthew G Knepley     Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
2354d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
2355d32f9abdSBarry Smith 
2356d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
2357d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
2358d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
2359d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
2360d32f9abdSBarry Smith 
2361c8a0d604SMatthew G Knepley $     For the Schur complement preconditioner if J = ( A00 A01 )
2362c8a0d604SMatthew G Knepley $                                                    ( A10 A11 )
2363c8a0d604SMatthew G Knepley $     the preconditioner using full factorization is
236413b05affSJed Brown $              ( I   -ksp(A00) A01 ) ( inv(A00)     0  ) (     I          0  )
2365c8a0d604SMatthew G Knepley $              ( 0         I       ) (   0      ksp(S) ) ( -A10 ksp(A00)  I  )
236613b05affSJed Brown      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_.  S is the Schur complement
236713b05affSJed Brown $              S = A11 - A10 ksp(A00) A01
236813b05affSJed Brown      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
2369c8a0d604SMatthew G Knepley      in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0,
2370c8a0d604SMatthew G Knepley      it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
2371a6a584a2SBarry Smith      A11 is used to construct a preconditioner for S, use PCFieldSplitSetSchurPre() for all the possible ways to construct the preconditioner for S.
23721d71051fSDmitry Karpeev 
2373a7476a74SDmitry Karpeev      The factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above,
23745668aaf4SBarry Smith      diag gives
2375c8a0d604SMatthew G Knepley $              ( inv(A00)     0   )
2376c8a0d604SMatthew G Knepley $              (   0      -ksp(S) )
23775668aaf4SBarry Smith      note that slightly counter intuitively there is a negative in front of the ksp(S) so that the preconditioner is positive definite. The lower factorization is the inverse of
2378c8a0d604SMatthew G Knepley $              (  A00   0 )
2379c8a0d604SMatthew G Knepley $              (  A10   S )
2380c8a0d604SMatthew G Knepley      where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of
2381c8a0d604SMatthew G Knepley $              ( A00 A01 )
2382c8a0d604SMatthew G Knepley $              (  0   S  )
2383c8a0d604SMatthew G Knepley      where again the inverses of A00 and S are applied using KSPs.
2384e69d4d44SBarry Smith 
2385edf189efSBarry Smith      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
2386edf189efSBarry Smith      is used automatically for a second block.
2387edf189efSBarry Smith 
2388ff218e97SBarry Smith      The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1.
2389ff218e97SBarry Smith      Generally it should be used with the AIJ format.
2390ff218e97SBarry Smith 
2391ff218e97SBarry Smith      The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see,
2392ff218e97SBarry Smith      for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT
2393ff218e97SBarry Smith      inside a smoother resulting in "Distributive Smoothers".
23940716a85fSBarry Smith 
2395a541d17aSBarry Smith    Concepts: physics based preconditioners, block preconditioners
23960971522cSBarry Smith 
23973bc631e0SBarry Smith    There is a nice discussion of block preconditioners in
23983bc631e0SBarry Smith 
23993bc631e0SBarry Smith [El08] A taxonomy and comparison of parallel block multi-level preconditioners for the incompressible Navier-Stokes equations
24003bc631e0SBarry Smith        Howard Elman, V.E. Howle, John Shadid, Robert Shuttleworth, Ray Tuminaro, Journal of Computational Physics 227 (2008) 1790--1808
24013bc631e0SBarry Smith        http://chess.cs.umd.edu/~elman/papers/tax.pdf
24023bc631e0SBarry Smith 
2403a6a584a2SBarry Smith    The Constrained Pressure Preconditioner (CPR) does not appear to be currently implementable directly with PCFIELDSPLIT. CPR solves first the Schur complemented pressure equation, updates the
2404a6a584a2SBarry Smith    residual on all variables and then applies a simple ILU like preconditioner on all the variables. So it is very much like the full Schur complement with selfp representing the Schur complement but instead
2405a6a584a2SBarry Smith    of backsolving for the saturations in the last step it solves a full coupled (ILU) system for updates to all the variables.
2406a6a584a2SBarry Smith 
24077e8cb189SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
24081d71051fSDmitry Karpeev            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSetSchurPre(),
24091d71051fSDmitry Karpeev 	   MatSchurComplementSetAinvType()
24100971522cSBarry Smith M*/
24110971522cSBarry Smith 
24120971522cSBarry Smith #undef __FUNCT__
24130971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
24148cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_FieldSplit(PC pc)
24150971522cSBarry Smith {
24160971522cSBarry Smith   PetscErrorCode ierr;
24170971522cSBarry Smith   PC_FieldSplit  *jac;
24180971522cSBarry Smith 
24190971522cSBarry Smith   PetscFunctionBegin;
2420b00a9115SJed Brown   ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr);
24212fa5cd67SKarl Rupp 
24220971522cSBarry Smith   jac->bs                 = -1;
24230971522cSBarry Smith   jac->nsplits            = 0;
24243e197d65SBarry Smith   jac->type               = PC_COMPOSITE_MULTIPLICATIVE;
2425e6cab6aaSJed Brown   jac->schurpre           = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
2426c9c6ffaaSJed Brown   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL;
2427fbe7908bSJed Brown   jac->dm_splits          = PETSC_TRUE;
242851f519a2SBarry Smith 
24290971522cSBarry Smith   pc->data = (void*)jac;
24300971522cSBarry Smith 
24310971522cSBarry Smith   pc->ops->apply           = PCApply_FieldSplit;
2432421e10b8SBarry Smith   pc->ops->applytranspose  = PCApplyTranspose_FieldSplit;
24330971522cSBarry Smith   pc->ops->setup           = PCSetUp_FieldSplit;
2434574deadeSBarry Smith   pc->ops->reset           = PCReset_FieldSplit;
24350971522cSBarry Smith   pc->ops->destroy         = PCDestroy_FieldSplit;
24360971522cSBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_FieldSplit;
24370971522cSBarry Smith   pc->ops->view            = PCView_FieldSplit;
24380971522cSBarry Smith   pc->ops->applyrichardson = 0;
24390971522cSBarry Smith 
2440bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
2441bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
2442bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
2443bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
2444bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
24456dbb499eSCian Wilson   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitRestrictIS_C",PCFieldSplitRestrictIS_FieldSplit);CHKERRQ(ierr);
24460971522cSBarry Smith   PetscFunctionReturn(0);
24470971522cSBarry Smith }
24480971522cSBarry Smith 
24490971522cSBarry Smith 
2450a541d17aSBarry Smith 
2451