xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 2bf68e3e0f2a61f71e7c65bee250bfa1c8ce0cdb)
1dba47a55SKris Buschelman 
2f5f0d762SBarry Smith 
3af0996ceSBarry Smith #include <petsc/private/pcimpl.h>     /*I "petscpc.h" I*/
4a80b646eSBarry Smith #include <petsc/private/kspimpl.h>    /*  This is needed to provide the appropriate PETSC_EXTERN for KSP_Solve_FS ....*/
51e25c274SJed Brown #include <petscdm.h>
60971522cSBarry Smith 
72e71c61dSMatthew G. Knepley const char *const PCFieldSplitSchurPreTypes[] = {"SELF","SELFP","A11","USER","FULL","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0};
8c9c6ffaaSJed Brown const char *const PCFieldSplitSchurFactTypes[] = {"DIAG","LOWER","UPPER","FULL","PCFieldSplitSchurFactType","PC_FIELDSPLIT_SCHUR_FACT_",0};
9c5d2311dSJed Brown 
109df09d43SBarry 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;
119df09d43SBarry Smith 
120971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
130971522cSBarry Smith struct _PC_FieldSplitLink {
1469a612a9SBarry Smith   KSP               ksp;
15443836d0SMatthew G Knepley   Vec               x,y,z;
16db4c96c1SJed Brown   char              *splitname;
170971522cSBarry Smith   PetscInt          nfields;
185d4c12cdSJungho Lee   PetscInt          *fields,*fields_col;
191b9fc7fcSBarry Smith   VecScatter        sctx;
20f5f0d762SBarry Smith   IS                is,is_col;
2151f519a2SBarry Smith   PC_FieldSplitLink next,previous;
229df09d43SBarry Smith   PetscLogEvent     event;
230971522cSBarry Smith };
240971522cSBarry Smith 
250971522cSBarry Smith typedef struct {
2668dd23aaSBarry Smith   PCCompositeType type;
27ace3abfcSBarry Smith   PetscBool       defaultsplit;                    /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */
28ace3abfcSBarry Smith   PetscBool       splitdefined;                    /* Flag is set after the splits have been defined, to prevent more splits from being added */
2930ad9308SMatthew Knepley   PetscInt        bs;                              /* Block size for IS and Mat structures */
3030ad9308SMatthew Knepley   PetscInt        nsplits;                         /* Number of field divisions defined */
3179416396SBarry Smith   Vec             *x,*y,w1,w2;
32519d70e2SJed Brown   Mat             *mat;                            /* The diagonal block for each split */
33519d70e2SJed Brown   Mat             *pmat;                           /* The preconditioning diagonal block for each split */
3430ad9308SMatthew Knepley   Mat             *Afield;                         /* The rows of the matrix associated with each split */
35ace3abfcSBarry Smith   PetscBool       issetup;
362fa5cd67SKarl Rupp 
3730ad9308SMatthew Knepley   /* Only used when Schur complement preconditioning is used */
3830ad9308SMatthew Knepley   Mat                       B;                     /* The (0,1) block */
3930ad9308SMatthew Knepley   Mat                       C;                     /* The (1,0) block */
40443836d0SMatthew G Knepley   Mat                       schur;                 /* The Schur complement S = A11 - A10 A00^{-1} A01, the KSP here, kspinner, is H_1 in [El08] */
41a7476a74SDmitry Karpeev   Mat                       schurp;                /* Assembled approximation to S built by MatSchurComplement to be used as a preconditioning matrix when solving with S */
42084e4875SJed Brown   Mat                       schur_user;            /* User-provided preconditioning matrix for the Schur complement */
43084e4875SJed Brown   PCFieldSplitSchurPreType  schurpre;              /* Determines which preconditioning matrix is used for the Schur complement */
44c9c6ffaaSJed Brown   PCFieldSplitSchurFactType schurfactorization;
4530ad9308SMatthew Knepley   KSP                       kspschur;              /* The solver for S */
46443836d0SMatthew G Knepley   KSP                       kspupper;              /* The solver for A in the upper diagonal part of the factorization (H_2 in [El08]) */
4797bbdb24SBarry Smith   PC_FieldSplitLink         head;
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   /*
327f5f0d762SBarry Smith    Kinda messy, but at least this now uses DMCreateFieldDecomposition().
3287287d2a3SDmitry Karpeev    Should probably be rewritten.
3297287d2a3SDmitry Karpeev    */
330f5f0d762SBarry Smith   if (!ilink) {
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);
4106ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
4116ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr);
412fcfd50ebSBarry Smith         ierr = ISDestroy(&zerodiags);CHKERRQ(ierr);
413fcfd50ebSBarry Smith         ierr = ISDestroy(&rest);CHKERRQ(ierr);
4143a062f41SBarry Smith       } else if (coupling) {
4153a062f41SBarry Smith         IS       coupling,rest;
4163a062f41SBarry Smith         PetscInt nmin,nmax;
4173a062f41SBarry Smith 
4183a062f41SBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
4193a062f41SBarry Smith         ierr = MatFindOffBlockDiagonalEntries(pc->mat,&coupling);CHKERRQ(ierr);
4206bea0878SBarry Smith         ierr = ISCreateStride(PetscObjectComm((PetscObject)pc->mat),nmax-nmin,nmin,1,&rest);CHKERRQ(ierr);
421e52d2c62SBarry Smith         ierr = ISSetIdentity(rest);CHKERRQ(ierr);
422d020c841SBarry Smith         ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
423d020c841SBarry Smith         ierr = PCFieldSplitSetIS(pc,"1",coupling);CHKERRQ(ierr);
4243a062f41SBarry Smith         ierr = ISDestroy(&coupling);CHKERRQ(ierr);
4253a062f41SBarry Smith         ierr = ISDestroy(&rest);CHKERRQ(ierr);
4266ce1633cSBarry Smith       } else {
427c5929fdfSBarry Smith         ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&fieldsplit_default,NULL);CHKERRQ(ierr);
4287287d2a3SDmitry Karpeev         if (!fieldsplit_default) {
429d32f9abdSBarry Smith           /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
430d32f9abdSBarry Smith            then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
4316c924f48SJed Brown           ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
4326c924f48SJed Brown           if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
433d32f9abdSBarry Smith         }
4346dbb499eSCian Wilson         if ((fieldsplit_default || !jac->splitdefined) && !jac->isrestrict) {
435d32f9abdSBarry Smith           ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
436db4c96c1SJed Brown           for (i=0; i<jac->bs; i++) {
4376c924f48SJed Brown             char splitname[8];
4388caf3d72SBarry Smith             ierr = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr);
4395d4c12cdSJungho Lee             ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr);
44079416396SBarry Smith           }
4415d4c12cdSJungho Lee           jac->defaultsplit = PETSC_TRUE;
442521d7252SBarry Smith         }
44366ffff09SJed Brown       }
4446ce1633cSBarry Smith     }
445edf189efSBarry Smith   } else if (jac->nsplits == 1) {
446edf189efSBarry Smith     if (ilink->is) {
447edf189efSBarry Smith       IS       is2;
448edf189efSBarry Smith       PetscInt nmin,nmax;
449edf189efSBarry Smith 
450edf189efSBarry Smith       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
451edf189efSBarry Smith       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
452db4c96c1SJed Brown       ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr);
453fcfd50ebSBarry Smith       ierr = ISDestroy(&is2);CHKERRQ(ierr);
45482f516ccSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
455edf189efSBarry Smith   }
456d0af7cd3SBarry Smith 
45782f516ccSBarry Smith   if (jac->nsplits < 2) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits);
45869a612a9SBarry Smith   PetscFunctionReturn(0);
45969a612a9SBarry Smith }
46069a612a9SBarry Smith 
461c5929fdfSBarry Smith PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(PetscOptions,const char pre[], const char name[], char *value[], PetscBool *flg);
462514bf10dSMatthew G Knepley 
46369a612a9SBarry Smith #undef __FUNCT__
46469a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
46569a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
46669a612a9SBarry Smith {
46769a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
46869a612a9SBarry Smith   PetscErrorCode    ierr;
4695a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
4702c9966d7SBarry Smith   PetscInt          i,nsplit;
4715f4ab4e1SJungho Lee   PetscBool         sorted, sorted_col;
47269a612a9SBarry Smith 
47369a612a9SBarry Smith   PetscFunctionBegin;
47469a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
47597bbdb24SBarry Smith   nsplit = jac->nsplits;
4765a9f2f41SSatish Balay   ilink  = jac->head;
47797bbdb24SBarry Smith 
47897bbdb24SBarry Smith   /* get the matrices for each split */
479704ba839SBarry Smith   if (!jac->issetup) {
4801b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
48197bbdb24SBarry Smith 
482704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
483704ba839SBarry Smith 
4845d4c12cdSJungho Lee     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
4852c9966d7SBarry Smith     if (jac->defaultsplit || !ilink->is) {
4862c9966d7SBarry Smith       if (jac->bs <= 0) jac->bs = nsplit;
4872c9966d7SBarry Smith     }
48851f519a2SBarry Smith     bs     = jac->bs;
48997bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
4901b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
4911b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
4921b9fc7fcSBarry Smith       if (jac->defaultsplit) {
493ce94432eSBarry Smith         ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
4945f4ab4e1SJungho Lee         ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr);
495704ba839SBarry Smith       } else if (!ilink->is) {
496ccb205f8SBarry Smith         if (ilink->nfields > 1) {
4975f4ab4e1SJungho Lee           PetscInt *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col;
498785e854fSJed Brown           ierr = PetscMalloc1(ilink->nfields*nslots,&ii);CHKERRQ(ierr);
499785e854fSJed Brown           ierr = PetscMalloc1(ilink->nfields*nslots,&jj);CHKERRQ(ierr);
5001b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
5011b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
5021b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
5035f4ab4e1SJungho Lee               jj[nfields*j + k] = rstart + bs*j + fields_col[k];
50497bbdb24SBarry Smith             }
50597bbdb24SBarry Smith           }
506ce94432eSBarry Smith           ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr);
507ce94432eSBarry Smith           ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr);
50890e68f20SPatrick Farrell           ierr = ISSetBlockSize(ilink->is, nfields);CHKERRQ(ierr);
50990e68f20SPatrick Farrell           ierr = ISSetBlockSize(ilink->is_col, nfields);CHKERRQ(ierr);
510ccb205f8SBarry Smith         } else {
511ce94432eSBarry Smith           ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
512ce94432eSBarry Smith           ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr);
513ccb205f8SBarry Smith        }
5143e197d65SBarry Smith       }
515704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
5165f4ab4e1SJungho Lee       if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); }
5175f4ab4e1SJungho Lee       if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
518704ba839SBarry Smith       ilink = ilink->next;
5191b9fc7fcSBarry Smith     }
5201b9fc7fcSBarry Smith   }
5211b9fc7fcSBarry Smith 
522704ba839SBarry Smith   ilink = jac->head;
52397bbdb24SBarry Smith   if (!jac->pmat) {
524bdddcaaaSMatthew G Knepley     Vec xtmp;
525bdddcaaaSMatthew G Knepley 
5262a7a6963SBarry Smith     ierr = MatCreateVecs(pc->pmat,&xtmp,NULL);CHKERRQ(ierr);
527785e854fSJed Brown     ierr = PetscMalloc1(nsplit,&jac->pmat);CHKERRQ(ierr);
528dcca6d9dSJed Brown     ierr = PetscMalloc2(nsplit,&jac->x,nsplit,&jac->y);CHKERRQ(ierr);
529cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
530bdddcaaaSMatthew G Knepley       MatNullSpace sp;
531bdddcaaaSMatthew G Knepley 
532a3df900dSMatthew G Knepley       /* Check for preconditioning matrix attached to IS */
533a3df900dSMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &jac->pmat[i]);CHKERRQ(ierr);
534a3df900dSMatthew G Knepley       if (jac->pmat[i]) {
535a3df900dSMatthew G Knepley         ierr = PetscObjectReference((PetscObject) jac->pmat[i]);CHKERRQ(ierr);
536a3df900dSMatthew G Knepley         if (jac->type == PC_COMPOSITE_SCHUR) {
537a3df900dSMatthew G Knepley           jac->schur_user = jac->pmat[i];
5382fa5cd67SKarl Rupp 
539a3df900dSMatthew G Knepley           ierr = PetscObjectReference((PetscObject) jac->schur_user);CHKERRQ(ierr);
540a3df900dSMatthew G Knepley         }
541a3df900dSMatthew G Knepley       } else {
5423a062f41SBarry Smith         const char *prefix;
5435f4ab4e1SJungho Lee         ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
5443a062f41SBarry Smith         ierr = KSPGetOptionsPrefix(ilink->ksp,&prefix);CHKERRQ(ierr);
5453a062f41SBarry Smith         ierr = MatSetOptionsPrefix(jac->pmat[i],prefix);CHKERRQ(ierr);
5463a062f41SBarry Smith         ierr = MatViewFromOptions(jac->pmat[i],NULL,"-mat_view");CHKERRQ(ierr);
547a3df900dSMatthew G Knepley       }
548bdddcaaaSMatthew G Knepley       /* create work vectors for each split */
5492a7a6963SBarry Smith       ierr     = MatCreateVecs(jac->pmat[i],&jac->x[i],&jac->y[i]);CHKERRQ(ierr);
5500298fd71SBarry Smith       ilink->x = jac->x[i]; ilink->y = jac->y[i]; ilink->z = NULL;
551bdddcaaaSMatthew G Knepley       /* compute scatter contexts needed by multiplicative versions and non-default splits */
5520298fd71SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],NULL,&ilink->sctx);CHKERRQ(ierr);
553ed1f0337SMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "nearnullspace", (PetscObject*) &sp);CHKERRQ(ierr);
554ed1f0337SMatthew G Knepley       if (sp) {
555ed1f0337SMatthew G Knepley         ierr = MatSetNearNullSpace(jac->pmat[i], sp);CHKERRQ(ierr);
556ed1f0337SMatthew G Knepley       }
557704ba839SBarry Smith       ilink = ilink->next;
558cf502942SBarry Smith     }
559bdddcaaaSMatthew G Knepley     ierr = VecDestroy(&xtmp);CHKERRQ(ierr);
56097bbdb24SBarry Smith   } else {
561cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
562a3df900dSMatthew G Knepley       Mat pmat;
563a3df900dSMatthew G Knepley 
564a3df900dSMatthew G Knepley       /* Check for preconditioning matrix attached to IS */
565a3df900dSMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &pmat);CHKERRQ(ierr);
566a3df900dSMatthew G Knepley       if (!pmat) {
5675f4ab4e1SJungho Lee         ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
568a3df900dSMatthew G Knepley       }
569704ba839SBarry Smith       ilink = ilink->next;
570cf502942SBarry Smith     }
57197bbdb24SBarry Smith   }
5724e39094bSDmitry Karpeev   if (jac->diag_use_amat) {
573519d70e2SJed Brown     ilink = jac->head;
574519d70e2SJed Brown     if (!jac->mat) {
575785e854fSJed Brown       ierr = PetscMalloc1(nsplit,&jac->mat);CHKERRQ(ierr);
576519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
5775f4ab4e1SJungho Lee         ierr  = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
578519d70e2SJed Brown         ilink = ilink->next;
579519d70e2SJed Brown       }
580519d70e2SJed Brown     } else {
581519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
5825f4ab4e1SJungho Lee         if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);}
583519d70e2SJed Brown         ilink = ilink->next;
584519d70e2SJed Brown       }
585519d70e2SJed Brown     }
586519d70e2SJed Brown   } else {
587519d70e2SJed Brown     jac->mat = jac->pmat;
588519d70e2SJed Brown   }
58997bbdb24SBarry Smith 
59053935eafSBarry Smith   /* Check for null space attached to IS */
59153935eafSBarry Smith   ilink = jac->head;
59253935eafSBarry Smith   for (i=0; i<nsplit; i++) {
59353935eafSBarry Smith     MatNullSpace sp;
59453935eafSBarry Smith 
59553935eafSBarry Smith     ierr = PetscObjectQuery((PetscObject) ilink->is, "nullspace", (PetscObject*) &sp);CHKERRQ(ierr);
59653935eafSBarry Smith     if (sp) {
59753935eafSBarry Smith       ierr = MatSetNullSpace(jac->mat[i], sp);CHKERRQ(ierr);
59853935eafSBarry Smith     }
59953935eafSBarry Smith     ilink = ilink->next;
60053935eafSBarry Smith   }
60153935eafSBarry Smith 
6026c8605c2SJed Brown   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
60368dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
6044e39094bSDmitry Karpeev     /* FIXME: Can/should we reuse jac->mat whenever (jac->diag_use_amat) is true? */
60568dd23aaSBarry Smith     ilink = jac->head;
606e52d2c62SBarry Smith     if (nsplit == 2 && jac->type == PC_COMPOSITE_MULTIPLICATIVE) {
607e52d2c62SBarry 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 */
608e52d2c62SBarry Smith       if (!jac->Afield) {
609e52d2c62SBarry Smith         ierr = PetscCalloc1(nsplit,&jac->Afield);CHKERRQ(ierr);
61080c96bb1SFande Kong         if (jac->offdiag_use_amat) {
611e52d2c62SBarry Smith           ierr  = MatGetSubMatrix(pc->mat,ilink->next->is,ilink->is,MAT_INITIAL_MATRIX,&jac->Afield[1]);CHKERRQ(ierr);
612e52d2c62SBarry Smith         } else {
61380c96bb1SFande Kong           ierr  = MatGetSubMatrix(pc->pmat,ilink->next->is,ilink->is,MAT_INITIAL_MATRIX,&jac->Afield[1]);CHKERRQ(ierr);
61480c96bb1SFande Kong         }
61580c96bb1SFande Kong       } else {
61680c96bb1SFande Kong         if (jac->offdiag_use_amat) {
617e52d2c62SBarry Smith           ierr  = MatGetSubMatrix(pc->mat,ilink->next->is,ilink->is,MAT_REUSE_MATRIX,&jac->Afield[1]);CHKERRQ(ierr);
61880c96bb1SFande Kong         } else {
61980c96bb1SFande Kong           ierr  = MatGetSubMatrix(pc->pmat,ilink->next->is,ilink->is,MAT_REUSE_MATRIX,&jac->Afield[1]);CHKERRQ(ierr);
62080c96bb1SFande Kong         }
621e52d2c62SBarry Smith       }
622e52d2c62SBarry Smith     } else {
62368dd23aaSBarry Smith       if (!jac->Afield) {
624785e854fSJed Brown         ierr = PetscMalloc1(nsplit,&jac->Afield);CHKERRQ(ierr);
62568dd23aaSBarry Smith         for (i=0; i<nsplit; i++) {
62680c96bb1SFande Kong           if (jac->offdiag_use_amat) {
6270298fd71SBarry Smith             ierr  = MatGetSubMatrix(pc->mat,ilink->is,NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
62880c96bb1SFande Kong           } else {
62980c96bb1SFande Kong             ierr  = MatGetSubMatrix(pc->pmat,ilink->is,NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
63080c96bb1SFande Kong           }
63168dd23aaSBarry Smith           ilink = ilink->next;
63268dd23aaSBarry Smith         }
63368dd23aaSBarry Smith       } else {
63468dd23aaSBarry Smith         for (i=0; i<nsplit; i++) {
63580c96bb1SFande Kong           if (jac->offdiag_use_amat) {
6360298fd71SBarry Smith             ierr  = MatGetSubMatrix(pc->mat,ilink->is,NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
63780c96bb1SFande Kong           } else {
63880c96bb1SFande Kong             ierr  = MatGetSubMatrix(pc->pmat,ilink->is,NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
63980c96bb1SFande Kong           }
64068dd23aaSBarry Smith           ilink = ilink->next;
64168dd23aaSBarry Smith         }
64268dd23aaSBarry Smith       }
64368dd23aaSBarry Smith     }
644e52d2c62SBarry Smith   }
64568dd23aaSBarry Smith 
6463b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
6473b224e63SBarry Smith     IS          ccis;
6484aa3045dSJed Brown     PetscInt    rstart,rend;
649093c86ffSJed Brown     char        lscname[256];
650093c86ffSJed Brown     PetscObject LSC_L;
651ce94432eSBarry Smith 
652ce94432eSBarry Smith     if (nsplit != 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
65368dd23aaSBarry Smith 
654e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
655e6cab6aaSJed Brown     ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
656e6cab6aaSJed Brown 
6573b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
6583b224e63SBarry Smith     if (jac->schur) {
6590298fd71SBarry Smith       KSP kspA = jac->head->ksp, kspInner = NULL, kspUpper = jac->kspupper;
660443836d0SMatthew G Knepley 
661fb3147dbSMatthew G Knepley       ierr  = MatSchurComplementGetKSP(jac->schur, &kspInner);CHKERRQ(ierr);
6623b224e63SBarry Smith       ilink = jac->head;
66349bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
6644e39094bSDmitry Karpeev       if (jac->offdiag_use_amat) {
6654aa3045dSJed Brown 	ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
666c84da90fSDmitry Karpeev       } else {
667df06c9e3SDmitry Karpeev 	ierr  = MatGetSubMatrix(pc->pmat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
668c84da90fSDmitry Karpeev       }
669fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
6703b224e63SBarry Smith       ilink = ilink->next;
67149bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
6724e39094bSDmitry Karpeev       if (jac->offdiag_use_amat) {
6734aa3045dSJed Brown 	ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
674c84da90fSDmitry Karpeev       } else {
675df06c9e3SDmitry Karpeev 	ierr  = MatGetSubMatrix(pc->pmat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
676c84da90fSDmitry Karpeev       }
677fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
678aa6c7ce3SBarry Smith       ierr  = MatSchurComplementUpdateSubMatrices(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);CHKERRQ(ierr);
679a7476a74SDmitry Karpeev       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELFP) {
680a7476a74SDmitry Karpeev 	ierr = MatDestroy(&jac->schurp);CHKERRQ(ierr);
681a7476a74SDmitry Karpeev 	ierr = MatSchurComplementGetPmat(jac->schur,MAT_INITIAL_MATRIX,&jac->schurp);CHKERRQ(ierr);
682a7476a74SDmitry Karpeev       }
683443836d0SMatthew G Knepley       if (kspA != kspInner) {
68423ee1639SBarry Smith         ierr = KSPSetOperators(kspA,jac->mat[0],jac->pmat[0]);CHKERRQ(ierr);
685443836d0SMatthew G Knepley       }
686443836d0SMatthew G Knepley       if (kspUpper != kspA) {
68723ee1639SBarry Smith         ierr = KSPSetOperators(kspUpper,jac->mat[0],jac->pmat[0]);CHKERRQ(ierr);
688443836d0SMatthew G Knepley       }
68923ee1639SBarry Smith       ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac));CHKERRQ(ierr);
6903b224e63SBarry Smith     } else {
691bafc1b83SMatthew G Knepley       const char   *Dprefix;
692470b340bSDmitry Karpeev       char         schurprefix[256], schurmatprefix[256];
693514bf10dSMatthew G Knepley       char         schurtestoption[256];
694bdddcaaaSMatthew G Knepley       MatNullSpace sp;
695514bf10dSMatthew G Knepley       PetscBool    flg;
6963b224e63SBarry Smith 
697a04f6461SBarry Smith       /* extract the A01 and A10 matrices */
6983b224e63SBarry Smith       ilink = jac->head;
69949bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
7004e39094bSDmitry Karpeev       if (jac->offdiag_use_amat) {
7014aa3045dSJed Brown 	ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
702c84da90fSDmitry Karpeev       } else {
703df06c9e3SDmitry Karpeev 	ierr  = MatGetSubMatrix(pc->pmat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
704c84da90fSDmitry Karpeev       }
705fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
7063b224e63SBarry Smith       ilink = ilink->next;
70749bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
7084e39094bSDmitry Karpeev       if (jac->offdiag_use_amat) {
7094aa3045dSJed Brown 	ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
710c84da90fSDmitry Karpeev       } else {
711df06c9e3SDmitry Karpeev 	ierr  = MatGetSubMatrix(pc->pmat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
712c84da90fSDmitry Karpeev       }
713fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
71420252d06SBarry Smith 
715f5236f50SJed Brown       /* Use mat[0] (diagonal block of Amat) preconditioned by pmat[0] to define Schur complement */
71620252d06SBarry Smith       ierr = MatCreate(((PetscObject)jac->mat[0])->comm,&jac->schur);CHKERRQ(ierr);
71720252d06SBarry Smith       ierr = MatSetType(jac->schur,MATSCHURCOMPLEMENT);CHKERRQ(ierr);
718bee83525SDmitry Karpeev       ierr = MatSchurComplementSetSubMatrices(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);CHKERRQ(ierr);
719470b340bSDmitry Karpeev       ierr = PetscSNPrintf(schurmatprefix, sizeof(schurmatprefix), "%sfieldsplit_%s_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr);
720470b340bSDmitry 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? */
721470b340bSDmitry Karpeev       ierr = MatSetOptionsPrefix(jac->schur,schurmatprefix);CHKERRQ(ierr);
722470b340bSDmitry Karpeev       ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
723895c21f2SBarry Smith       ierr = MatGetNullSpace(jac->mat[1], &sp);CHKERRQ(ierr);
72420252d06SBarry Smith       if (sp) {
72520252d06SBarry Smith         ierr = MatSetNullSpace(jac->schur, sp);CHKERRQ(ierr);
72620252d06SBarry Smith       }
72720252d06SBarry Smith 
72820252d06SBarry Smith       ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_inner_", ilink->splitname);CHKERRQ(ierr);
729c5929fdfSBarry Smith       ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->options,((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);CHKERRQ(ierr);
730514bf10dSMatthew G Knepley       if (flg) {
731514bf10dSMatthew G Knepley         DM  dmInner;
73221635b76SJed Brown         KSP kspInner;
73321635b76SJed Brown 
73421635b76SJed Brown         ierr = MatSchurComplementGetKSP(jac->schur, &kspInner);CHKERRQ(ierr);
73521635b76SJed Brown         ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_inner_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr);
73621635b76SJed Brown         /* Indent this deeper to emphasize the "inner" nature of this solver. */
73721635b76SJed Brown         ierr = PetscObjectIncrementTabLevel((PetscObject)kspInner, (PetscObject) pc, 2);CHKERRQ(ierr);
73821635b76SJed Brown         ierr = KSPSetOptionsPrefix(kspInner, schurprefix);CHKERRQ(ierr);
739514bf10dSMatthew G Knepley 
740514bf10dSMatthew G Knepley         /* Set DM for new solver */
741bafc1b83SMatthew G Knepley         ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr);
74221635b76SJed Brown         ierr = KSPSetDM(kspInner, dmInner);CHKERRQ(ierr);
74321635b76SJed Brown         ierr = KSPSetDMActive(kspInner, PETSC_FALSE);CHKERRQ(ierr);
744514bf10dSMatthew G Knepley       } else {
74521635b76SJed Brown          /* Use the outer solver for the inner solve, but revert the KSPPREONLY from PCFieldSplitSetFields_FieldSplit or
74621635b76SJed Brown           * PCFieldSplitSetIS_FieldSplit. We don't want KSPPREONLY because it makes the Schur complement inexact,
74721635b76SJed Brown           * preventing Schur complement reduction to be an accurate solve. Usually when an iterative solver is used for
74821635b76SJed 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
74921635b76SJed Brown           * GMRES the default. Note that it is also common to use PREONLY for S, in which case S may not be used
75021635b76SJed Brown           * directly, and the user is responsible for setting an inexact method for fieldsplit's A^{-1}. */
75121635b76SJed Brown         ierr = KSPSetType(jac->head->ksp,KSPGMRES);CHKERRQ(ierr);
752514bf10dSMatthew G Knepley         ierr = MatSchurComplementSetKSP(jac->schur,jac->head->ksp);CHKERRQ(ierr);
753bafc1b83SMatthew G Knepley       }
75423ee1639SBarry Smith       ierr = KSPSetOperators(jac->head->ksp,jac->mat[0],jac->pmat[0]);CHKERRQ(ierr);
7555a9f2f41SSatish Balay       ierr = KSPSetFromOptions(jac->head->ksp);CHKERRQ(ierr);
75697bbdb24SBarry Smith       ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
75797bbdb24SBarry Smith 
758443836d0SMatthew G Knepley       ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_upper_", ilink->splitname);CHKERRQ(ierr);
759c5929fdfSBarry Smith       ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->options,((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);CHKERRQ(ierr);
760443836d0SMatthew G Knepley       if (flg) {
761443836d0SMatthew G Knepley         DM dmInner;
762443836d0SMatthew G Knepley 
763443836d0SMatthew G Knepley         ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_upper_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr);
76482f516ccSBarry Smith         ierr = KSPCreate(PetscObjectComm((PetscObject)pc), &jac->kspupper);CHKERRQ(ierr);
765422a814eSBarry Smith         ierr = KSPSetErrorIfNotConverged(jac->kspupper,pc->erroriffailure);CHKERRQ(ierr);
766443836d0SMatthew G Knepley         ierr = KSPSetOptionsPrefix(jac->kspupper, schurprefix);CHKERRQ(ierr);
767443836d0SMatthew G Knepley         ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr);
768443836d0SMatthew G Knepley         ierr = KSPSetDM(jac->kspupper, dmInner);CHKERRQ(ierr);
769443836d0SMatthew G Knepley         ierr = KSPSetDMActive(jac->kspupper, PETSC_FALSE);CHKERRQ(ierr);
770443836d0SMatthew G Knepley         ierr = KSPSetFromOptions(jac->kspupper);CHKERRQ(ierr);
77123ee1639SBarry Smith         ierr = KSPSetOperators(jac->kspupper,jac->mat[0],jac->pmat[0]);CHKERRQ(ierr);
772443836d0SMatthew G Knepley         ierr = VecDuplicate(jac->head->x, &jac->head->z);CHKERRQ(ierr);
773443836d0SMatthew G Knepley       } else {
774443836d0SMatthew G Knepley         jac->kspupper = jac->head->ksp;
775443836d0SMatthew G Knepley         ierr          = PetscObjectReference((PetscObject) jac->head->ksp);CHKERRQ(ierr);
776443836d0SMatthew G Knepley       }
777443836d0SMatthew G Knepley 
778a7476a74SDmitry Karpeev       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELFP) {
779a7476a74SDmitry Karpeev 	ierr = MatSchurComplementGetPmat(jac->schur,MAT_INITIAL_MATRIX,&jac->schurp);CHKERRQ(ierr);
780a7476a74SDmitry Karpeev       }
781ce94432eSBarry Smith       ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&jac->kspschur);CHKERRQ(ierr);
782422a814eSBarry Smith       ierr = KSPSetErrorIfNotConverged(jac->kspschur,pc->erroriffailure);CHKERRQ(ierr);
78397bbdb24SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr);
78420252d06SBarry Smith       ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
78597bbdb24SBarry Smith       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
7867233a360SDmitry Karpeev         PC pcschur;
7877233a360SDmitry Karpeev         ierr = KSPGetPC(jac->kspschur,&pcschur);CHKERRQ(ierr);
7887233a360SDmitry Karpeev         ierr = PCSetType(pcschur,PCNONE);CHKERRQ(ierr);
78997bbdb24SBarry Smith         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
790e74569cdSMatthew G. Knepley       } else if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_FULL) {
791e74569cdSMatthew G. Knepley         ierr = MatSchurComplementComputeExplicitOperator(jac->schur, &jac->schur_user);CHKERRQ(ierr);
79297bbdb24SBarry Smith       }
79323ee1639SBarry Smith       ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac));CHKERRQ(ierr);
79497bbdb24SBarry Smith       ierr = KSPGetOptionsPrefix(jac->head->next->ksp, &Dprefix);CHKERRQ(ierr);
79597bbdb24SBarry Smith       ierr = KSPSetOptionsPrefix(jac->kspschur,         Dprefix);CHKERRQ(ierr);
796b20b4189SMatthew G. Knepley       /* propogate DM */
797b20b4189SMatthew G. Knepley       {
798b20b4189SMatthew G. Knepley         DM sdm;
799b20b4189SMatthew G. Knepley         ierr = KSPGetDM(jac->head->next->ksp, &sdm);CHKERRQ(ierr);
800b20b4189SMatthew G. Knepley         if (sdm) {
801b20b4189SMatthew G. Knepley           ierr = KSPSetDM(jac->kspschur, sdm);CHKERRQ(ierr);
802b20b4189SMatthew G. Knepley           ierr = KSPSetDMActive(jac->kspschur, PETSC_FALSE);CHKERRQ(ierr);
803b20b4189SMatthew G. Knepley         }
804b20b4189SMatthew G. Knepley       }
80597bbdb24SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
80697bbdb24SBarry Smith       /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */
80797bbdb24SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
80897bbdb24SBarry Smith     }
80997bbdb24SBarry Smith 
8105a9f2f41SSatish Balay     /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */
8118caf3d72SBarry Smith     ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_L",ilink->splitname);CHKERRQ(ierr);
812519d70e2SJed Brown     ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
8133b224e63SBarry Smith     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
814c1570756SJed Brown     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);}
8158caf3d72SBarry Smith     ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr);
81697bbdb24SBarry Smith     ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
8175a9f2f41SSatish Balay     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
8180971522cSBarry Smith     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);}
81997bbdb24SBarry Smith   } else {
82068bd789dSDmitry Karpeev     /* set up the individual splits' PCs */
8210971522cSBarry Smith     i     = 0;
8220971522cSBarry Smith     ilink = jac->head;
8230971522cSBarry Smith     while (ilink) {
82423ee1639SBarry Smith       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i]);CHKERRQ(ierr);
8250971522cSBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
8260971522cSBarry Smith       if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);}
8270971522cSBarry Smith       i++;
8280971522cSBarry Smith       ilink = ilink->next;
8290971522cSBarry Smith     }
8303b224e63SBarry Smith   }
8313b224e63SBarry Smith 
832c1570756SJed Brown   jac->suboptionsset = PETSC_TRUE;
8330971522cSBarry Smith   PetscFunctionReturn(0);
8340971522cSBarry Smith }
8350971522cSBarry Smith 
8365a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
837ca9f406cSSatish Balay   (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
838ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
8399df09d43SBarry Smith    PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL) ||\
8405a9f2f41SSatish Balay    KSPSolve(ilink->ksp,ilink->x,ilink->y) ||                               \
8419df09d43SBarry Smith    PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL) ||\
842ca9f406cSSatish Balay    VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) ||  \
843ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
84479416396SBarry Smith 
8450971522cSBarry Smith #undef __FUNCT__
8463b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur"
8473b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
8483b224e63SBarry Smith {
8493b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
8503b224e63SBarry Smith   PetscErrorCode    ierr;
8513b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
852443836d0SMatthew G Knepley   KSP               kspA   = ilinkA->ksp, kspLower = kspA, kspUpper = jac->kspupper;
8533b224e63SBarry Smith 
8543b224e63SBarry Smith   PetscFunctionBegin;
855c5d2311dSJed Brown   switch (jac->schurfactorization) {
856c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_DIAG:
857a04f6461SBarry Smith     /* [A00 0; 0 -S], positive definite, suitable for MINRES */
858c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
859c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
860c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8619df09d43SBarry Smith     ierr = PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
862443836d0SMatthew G Knepley     ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
8639df09d43SBarry Smith     ierr = PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
864c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
865c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8669df09d43SBarry Smith     ierr = PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
867c5d2311dSJed Brown     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
8689df09d43SBarry Smith     ierr = PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
869c5d2311dSJed Brown     ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr);
870c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
871c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
872c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
873c5d2311dSJed Brown     break;
874c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_LOWER:
875a04f6461SBarry Smith     /* [A00 0; A10 S], suitable for left preconditioning */
876c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
877c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8789df09d43SBarry Smith     ierr = PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
879443836d0SMatthew G Knepley     ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
8809df09d43SBarry Smith     ierr = PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
881c5d2311dSJed Brown     ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
882c5d2311dSJed Brown     ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr);
883c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
884c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
885c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8869df09d43SBarry Smith     ierr = PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
887c5d2311dSJed Brown     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
8889df09d43SBarry Smith     ierr = PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
889c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
890c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
891c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
892c5d2311dSJed Brown     break;
893c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_UPPER:
894a04f6461SBarry Smith     /* [A00 A01; 0 S], suitable for right preconditioning */
895c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
896c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8979df09d43SBarry Smith     ierr = PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
898c5d2311dSJed Brown     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
8999df09d43SBarry 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);
900c5d2311dSJed Brown     ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr);
901c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
902c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
903c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9049df09d43SBarry Smith     ierr = PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
905443836d0SMatthew G Knepley     ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
9069df09d43SBarry Smith     ierr = PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
907c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
908c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
909c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
910c5d2311dSJed Brown     break;
911c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_FULL:
912a04f6461SBarry 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 */
9133b224e63SBarry Smith     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9143b224e63SBarry Smith     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9159df09d43SBarry Smith     ierr = PetscLogEventBegin(KSP_Solve_FS_L,kspLower,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
916443836d0SMatthew G Knepley     ierr = KSPSolve(kspLower,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
9179df09d43SBarry Smith     ierr = PetscLogEventEnd(KSP_Solve_FS_L,kspLower,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
9183b224e63SBarry Smith     ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
9193b224e63SBarry Smith     ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
9203b224e63SBarry Smith     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9213b224e63SBarry Smith     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9223b224e63SBarry Smith 
9239df09d43SBarry Smith     ierr = PetscLogEventBegin(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
9243b224e63SBarry Smith     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
9259df09d43SBarry Smith     ierr = PetscLogEventEnd(KSP_Solve_FS_S,jac->kspschur,ilinkD->x,ilinkD->y,NULL);CHKERRQ(ierr);
9263b224e63SBarry Smith     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9273b224e63SBarry Smith     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9283b224e63SBarry Smith 
929443836d0SMatthew G Knepley     if (kspUpper == kspA) {
9303b224e63SBarry Smith       ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
9313b224e63SBarry Smith       ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
9329df09d43SBarry Smith       ierr = PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
933443836d0SMatthew G Knepley       ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
9349df09d43SBarry Smith       ierr = PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
935443836d0SMatthew G Knepley     } else {
9369df09d43SBarry Smith       ierr = PetscLogEventBegin(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
937443836d0SMatthew G Knepley       ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
9389df09d43SBarry Smith       ierr = PetscLogEventEnd(ilinkA->event,kspA,ilinkA->x,ilinkA->y,NULL);CHKERRQ(ierr);
939443836d0SMatthew G Knepley       ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
9409df09d43SBarry Smith       ierr = PetscLogEventBegin(KSP_Solve_FS_U,kspUpper,ilinkA->x,ilinkA->z,NULL);CHKERRQ(ierr);
941443836d0SMatthew G Knepley       ierr = KSPSolve(kspUpper,ilinkA->x,ilinkA->z);CHKERRQ(ierr);
9429df09d43SBarry Smith       ierr = PetscLogEventEnd(KSP_Solve_FS_U,kspUpper,ilinkA->x,ilinkA->z,NULL);CHKERRQ(ierr);
943443836d0SMatthew G Knepley       ierr = VecAXPY(ilinkA->y,-1.0,ilinkA->z);CHKERRQ(ierr);
944443836d0SMatthew G Knepley     }
9453b224e63SBarry Smith     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9463b224e63SBarry Smith     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
947c5d2311dSJed Brown   }
9483b224e63SBarry Smith   PetscFunctionReturn(0);
9493b224e63SBarry Smith }
9503b224e63SBarry Smith 
9513b224e63SBarry Smith #undef __FUNCT__
9520971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
9530971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
9540971522cSBarry Smith {
9550971522cSBarry Smith   PC_FieldSplit      *jac = (PC_FieldSplit*)pc->data;
9560971522cSBarry Smith   PetscErrorCode     ierr;
9575a9f2f41SSatish Balay   PC_FieldSplitLink  ilink = jac->head;
958939b8a20SBarry Smith   PetscInt           cnt,bs;
959568ed31eSHong Zhang   KSPConvergedReason reason;
9600971522cSBarry Smith 
9610971522cSBarry Smith   PetscFunctionBegin;
96279416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
9631b9fc7fcSBarry Smith     if (jac->defaultsplit) {
964939b8a20SBarry Smith       ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
965ce94432eSBarry 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);
966939b8a20SBarry Smith       ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr);
967ce94432eSBarry 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);
9680971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
9695a9f2f41SSatish Balay       while (ilink) {
9709df09d43SBarry Smith         ierr  = PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
9715a9f2f41SSatish Balay         ierr  = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
9729df09d43SBarry Smith         ierr  = PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
973568ed31eSHong Zhang         ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
974568ed31eSHong Zhang         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
975568ed31eSHong Zhang           pc->failedreason = PC_SUBPC_ERROR;
976568ed31eSHong Zhang         }
9775a9f2f41SSatish Balay         ilink = ilink->next;
9780971522cSBarry Smith       }
9790971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
9801b9fc7fcSBarry Smith     } else {
981efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
9825a9f2f41SSatish Balay       while (ilink) {
9835a9f2f41SSatish Balay         ierr  = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
984568ed31eSHong Zhang         ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
985568ed31eSHong Zhang         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
986568ed31eSHong Zhang           pc->failedreason = PC_SUBPC_ERROR;
987568ed31eSHong Zhang         }
9885a9f2f41SSatish Balay         ilink = ilink->next;
9891b9fc7fcSBarry Smith       }
9901b9fc7fcSBarry Smith     }
991e52d2c62SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE && jac->nsplits == 2) {
992e52d2c62SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
993e52d2c62SBarry Smith     /* solve on first block for first block variables */
994e52d2c62SBarry Smith     ierr = VecScatterBegin(ilink->sctx,x,ilink->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
995e52d2c62SBarry Smith     ierr = VecScatterEnd(ilink->sctx,x,ilink->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9969df09d43SBarry Smith     ierr = PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
997e52d2c62SBarry Smith     ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
9989df09d43SBarry Smith     ierr = PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
999568ed31eSHong Zhang     ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
1000568ed31eSHong Zhang     if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1001568ed31eSHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
1002568ed31eSHong Zhang     }
1003e52d2c62SBarry Smith     ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1004e52d2c62SBarry Smith     ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1005e52d2c62SBarry Smith 
1006e52d2c62SBarry Smith     /* compute the residual only onto second block variables using first block variables */
1007e52d2c62SBarry Smith     ierr = MatMult(jac->Afield[1],ilink->y,ilink->next->x);CHKERRQ(ierr);
1008e52d2c62SBarry Smith     ilink = ilink->next;
1009e52d2c62SBarry Smith     ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
1010e52d2c62SBarry Smith     ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1011e52d2c62SBarry Smith     ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1012e52d2c62SBarry Smith 
1013e52d2c62SBarry Smith     /* solve on second block variables */
10149df09d43SBarry Smith     ierr = PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
1015e52d2c62SBarry Smith     ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
10169df09d43SBarry Smith     ierr = PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
1017568ed31eSHong Zhang     ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
1018568ed31eSHong Zhang     if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1019568ed31eSHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
1020568ed31eSHong Zhang     }
1021e52d2c62SBarry Smith     ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1022e52d2c62SBarry Smith     ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
102316913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
102479416396SBarry Smith     if (!jac->w1) {
102579416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
102679416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
102779416396SBarry Smith     }
1028efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
10295a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
1030568ed31eSHong Zhang     ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
1031568ed31eSHong Zhang     if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1032568ed31eSHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
1033568ed31eSHong Zhang     }
10343e197d65SBarry Smith     cnt  = 1;
10355a9f2f41SSatish Balay     while (ilink->next) {
10365a9f2f41SSatish Balay       ilink = ilink->next;
10373e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
10383e197d65SBarry Smith       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
10393e197d65SBarry Smith       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
10403e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10413e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10429df09d43SBarry Smith       ierr = PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
10433e197d65SBarry Smith       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
10449df09d43SBarry Smith       ierr = PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
10452618eb8fSHong Zhang       ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
10462618eb8fSHong Zhang       if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
10472618eb8fSHong Zhang         pc->failedreason = PC_SUBPC_ERROR;
10482618eb8fSHong Zhang       }
10493e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10503e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10513e197d65SBarry Smith     }
105251f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
105311755939SBarry Smith       cnt -= 2;
105451f519a2SBarry Smith       while (ilink->previous) {
105551f519a2SBarry Smith         ilink = ilink->previous;
105611755939SBarry Smith         /* compute the residual only over the part of the vector needed */
105711755939SBarry Smith         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
105811755939SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
105911755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
106011755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10619df09d43SBarry Smith         ierr = PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
106211755939SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
10639df09d43SBarry Smith         ierr = PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
1064568ed31eSHong Zhang         ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
1065568ed31eSHong Zhang         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
1066568ed31eSHong Zhang           pc->failedreason = PC_SUBPC_ERROR;
1067568ed31eSHong Zhang         }
106811755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
106911755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
107051f519a2SBarry Smith       }
107111755939SBarry Smith     }
1072ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
10730971522cSBarry Smith   PetscFunctionReturn(0);
10740971522cSBarry Smith }
10750971522cSBarry Smith 
1076421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
1077ca9f406cSSatish Balay   (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
1078ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
10799df09d43SBarry Smith    PetscLogEventBegin(ilink->event,ilink->ksp,ilink->y,ilink->x,NULL) || \
1080421e10b8SBarry Smith    KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) ||                  \
10819df09d43SBarry Smith    PetscLogEventBegin(ilink->event,ilink->ksp,ilink->y,ilink->x,NULL) || \
1082ca9f406cSSatish Balay    VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
1083ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
1084421e10b8SBarry Smith 
1085421e10b8SBarry Smith #undef __FUNCT__
10868c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit"
1087421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
1088421e10b8SBarry Smith {
1089421e10b8SBarry Smith   PC_FieldSplit      *jac = (PC_FieldSplit*)pc->data;
1090421e10b8SBarry Smith   PetscErrorCode     ierr;
1091421e10b8SBarry Smith   PC_FieldSplitLink  ilink = jac->head;
1092939b8a20SBarry Smith   PetscInt           bs;
1093291dd214SHong Zhang   KSPConvergedReason reason;
1094421e10b8SBarry Smith 
1095421e10b8SBarry Smith   PetscFunctionBegin;
1096421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
1097421e10b8SBarry Smith     if (jac->defaultsplit) {
1098939b8a20SBarry Smith       ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
1099ce94432eSBarry 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);
1100939b8a20SBarry Smith       ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr);
1101ce94432eSBarry 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);
1102421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
1103421e10b8SBarry Smith       while (ilink) {
11049df09d43SBarry Smith         ierr = PetscLogEventBegin(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
1105421e10b8SBarry Smith         ierr  = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
11069df09d43SBarry Smith         ierr = PetscLogEventEnd(ilink->event,ilink->ksp,ilink->x,ilink->y,NULL);CHKERRQ(ierr);
11072618eb8fSHong Zhang         ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
11082618eb8fSHong Zhang         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
11092618eb8fSHong Zhang           pc->failedreason = PC_SUBPC_ERROR;
11102618eb8fSHong Zhang         }
1111421e10b8SBarry Smith         ilink = ilink->next;
1112421e10b8SBarry Smith       }
1113421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
1114421e10b8SBarry Smith     } else {
1115421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
1116421e10b8SBarry Smith       while (ilink) {
1117421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
11182618eb8fSHong Zhang         ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
11192618eb8fSHong Zhang         if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
11202618eb8fSHong Zhang           pc->failedreason = PC_SUBPC_ERROR;
11212618eb8fSHong Zhang         }
1122421e10b8SBarry Smith         ilink = ilink->next;
1123421e10b8SBarry Smith       }
1124421e10b8SBarry Smith     }
1125421e10b8SBarry Smith   } else {
1126421e10b8SBarry Smith     if (!jac->w1) {
1127421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
1128421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
1129421e10b8SBarry Smith     }
1130421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
1131421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
1132421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
11332618eb8fSHong Zhang       ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
11342618eb8fSHong Zhang       if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
11352618eb8fSHong Zhang         pc->failedreason = PC_SUBPC_ERROR;
11362618eb8fSHong Zhang       }
1137421e10b8SBarry Smith       while (ilink->next) {
1138421e10b8SBarry Smith         ilink = ilink->next;
11399989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
1140421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
1141421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
1142421e10b8SBarry Smith       }
1143421e10b8SBarry Smith       while (ilink->previous) {
1144421e10b8SBarry Smith         ilink = ilink->previous;
11459989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
1146421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
1147421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
1148421e10b8SBarry Smith       }
1149421e10b8SBarry Smith     } else {
1150421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
1151421e10b8SBarry Smith         ilink = ilink->next;
1152421e10b8SBarry Smith       }
1153421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
11542618eb8fSHong Zhang       ierr = KSPGetConvergedReason(ilink->ksp,&reason);CHKERRQ(ierr);
11552618eb8fSHong Zhang       if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
11562618eb8fSHong Zhang         pc->failedreason = PC_SUBPC_ERROR;
11572618eb8fSHong Zhang       }
1158421e10b8SBarry Smith       while (ilink->previous) {
1159421e10b8SBarry Smith         ilink = ilink->previous;
11609989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
1161421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
1162421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
1163421e10b8SBarry Smith       }
1164421e10b8SBarry Smith     }
1165421e10b8SBarry Smith   }
1166421e10b8SBarry Smith   PetscFunctionReturn(0);
1167421e10b8SBarry Smith }
1168421e10b8SBarry Smith 
11690971522cSBarry Smith #undef __FUNCT__
1170574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit"
1171574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc)
11720971522cSBarry Smith {
11730971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
11740971522cSBarry Smith   PetscErrorCode    ierr;
11755a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
11760971522cSBarry Smith 
11770971522cSBarry Smith   PetscFunctionBegin;
11785a9f2f41SSatish Balay   while (ilink) {
1179f5f0d762SBarry Smith     ierr  = KSPDestroy(&ilink->ksp);CHKERRQ(ierr);
1180fcfd50ebSBarry Smith     ierr  = VecDestroy(&ilink->x);CHKERRQ(ierr);
1181fcfd50ebSBarry Smith     ierr  = VecDestroy(&ilink->y);CHKERRQ(ierr);
1182443836d0SMatthew G Knepley     ierr  = VecDestroy(&ilink->z);CHKERRQ(ierr);
1183fcfd50ebSBarry Smith     ierr  = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr);
1184fcfd50ebSBarry Smith     ierr  = ISDestroy(&ilink->is);CHKERRQ(ierr);
1185b5787286SJed Brown     ierr  = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
1186f5f0d762SBarry Smith     ierr  = PetscFree(ilink->splitname);CHKERRQ(ierr);
1187f5f0d762SBarry Smith     ierr  = PetscFree(ilink->fields);CHKERRQ(ierr);
1188f5f0d762SBarry Smith     ierr  = PetscFree(ilink->fields_col);CHKERRQ(ierr);
11895a9f2f41SSatish Balay     next  = ilink->next;
1190f5f0d762SBarry Smith     ierr  = PetscFree(ilink);CHKERRQ(ierr);
11915a9f2f41SSatish Balay     ilink = next;
11920971522cSBarry Smith   }
1193f5f0d762SBarry Smith   jac->head = NULL;
119405b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
1195574deadeSBarry Smith   if (jac->mat && jac->mat != jac->pmat) {
1196574deadeSBarry Smith     ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);
1197574deadeSBarry Smith   } else if (jac->mat) {
11980298fd71SBarry Smith     jac->mat = NULL;
1199574deadeSBarry Smith   }
120097bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
120168dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
1202f5f0d762SBarry Smith   jac->nsplits = 0;
12036bf464f9SBarry Smith   ierr       = VecDestroy(&jac->w1);CHKERRQ(ierr);
12046bf464f9SBarry Smith   ierr       = VecDestroy(&jac->w2);CHKERRQ(ierr);
12056bf464f9SBarry Smith   ierr       = MatDestroy(&jac->schur);CHKERRQ(ierr);
1206a7476a74SDmitry Karpeev   ierr       = MatDestroy(&jac->schurp);CHKERRQ(ierr);
12076bf464f9SBarry Smith   ierr       = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
12086bf464f9SBarry Smith   ierr       = KSPDestroy(&jac->kspschur);CHKERRQ(ierr);
1209d78dad28SBarry Smith   ierr       = KSPDestroy(&jac->kspupper);CHKERRQ(ierr);
12106bf464f9SBarry Smith   ierr       = MatDestroy(&jac->B);CHKERRQ(ierr);
12116bf464f9SBarry Smith   ierr       = MatDestroy(&jac->C);CHKERRQ(ierr);
12126dbb499eSCian Wilson   jac->isrestrict = PETSC_FALSE;
1213574deadeSBarry Smith   PetscFunctionReturn(0);
1214574deadeSBarry Smith }
1215574deadeSBarry Smith 
1216574deadeSBarry Smith #undef __FUNCT__
1217574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
1218574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
1219574deadeSBarry Smith {
1220574deadeSBarry Smith   PetscErrorCode    ierr;
1221574deadeSBarry Smith 
1222574deadeSBarry Smith   PetscFunctionBegin;
1223574deadeSBarry Smith   ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr);
1224c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1225bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",NULL);CHKERRQ(ierr);
1226bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",NULL);CHKERRQ(ierr);
1227bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",NULL);CHKERRQ(ierr);
1228bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",NULL);CHKERRQ(ierr);
1229bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",NULL);CHKERRQ(ierr);
123029f8a81cSJed Brown   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",NULL);CHKERRQ(ierr);
123137a82bf0SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",NULL);CHKERRQ(ierr);
1232bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",NULL);CHKERRQ(ierr);
12336dbb499eSCian Wilson   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitRestrictIS_C",NULL);CHKERRQ(ierr);
12340971522cSBarry Smith   PetscFunctionReturn(0);
12350971522cSBarry Smith }
12360971522cSBarry Smith 
12370971522cSBarry Smith #undef __FUNCT__
12380971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
12394416b707SBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PetscOptionItems *PetscOptionsObject,PC pc)
12400971522cSBarry Smith {
12411b9fc7fcSBarry Smith   PetscErrorCode  ierr;
12426c924f48SJed Brown   PetscInt        bs;
1243bc59fbc5SBarry Smith   PetscBool       flg,stokes = PETSC_FALSE;
12449dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
12453b224e63SBarry Smith   PCCompositeType ctype;
12461b9fc7fcSBarry Smith 
12470971522cSBarry Smith   PetscFunctionBegin;
1248e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"FieldSplit options");CHKERRQ(ierr);
12494ab8060aSDmitry Karpeev   ierr = PetscOptionsBool("-pc_fieldsplit_dm_splits","Whether to use DMCreateFieldDecomposition() for splits","PCFieldSplitSetDMSplits",jac->dm_splits,&jac->dm_splits,NULL);CHKERRQ(ierr);
125051f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
125151f519a2SBarry Smith   if (flg) {
125251f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
125351f519a2SBarry Smith   }
12542686e3e9SMatthew G. Knepley   jac->diag_use_amat = pc->useAmat;
12552686e3e9SMatthew 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);
12562686e3e9SMatthew G. Knepley   jac->offdiag_use_amat = pc->useAmat;
12572686e3e9SMatthew 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);
125800216929SDmitry Karpeev   /* FIXME: No programmatic equivalent to the following. */
1259c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);CHKERRQ(ierr);
1260c0adfefeSBarry Smith   if (stokes) {
1261c0adfefeSBarry Smith     ierr          = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
1262c0adfefeSBarry Smith     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
1263c0adfefeSBarry Smith   }
1264c0adfefeSBarry Smith 
12653b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
12663b224e63SBarry Smith   if (flg) {
12673b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
12683b224e63SBarry Smith   }
1269c30613efSMatthew Knepley   /* Only setup fields once */
1270c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
1271d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
1272d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
12736c924f48SJed Brown     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
12746c924f48SJed Brown     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
1275d32f9abdSBarry Smith   }
1276c5d2311dSJed Brown   if (jac->type == PC_COMPOSITE_SCHUR) {
1277c5929fdfSBarry Smith     ierr = PetscOptionsGetEnum(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr);
1278c9c6ffaaSJed Brown     if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);}
12790298fd71SBarry 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);
128029f8a81cSJed 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);
1281c5d2311dSJed Brown   }
12821b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
12830971522cSBarry Smith   PetscFunctionReturn(0);
12840971522cSBarry Smith }
12850971522cSBarry Smith 
12860971522cSBarry Smith /*------------------------------------------------------------------------------------*/
12870971522cSBarry Smith 
12880971522cSBarry Smith #undef __FUNCT__
12890971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
12901e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
12910971522cSBarry Smith {
129297bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
12930971522cSBarry Smith   PetscErrorCode    ierr;
12945a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
129569a612a9SBarry Smith   char              prefix[128];
12965d4c12cdSJungho Lee   PetscInt          i;
12970971522cSBarry Smith 
12980971522cSBarry Smith   PetscFunctionBegin;
12996c924f48SJed Brown   if (jac->splitdefined) {
13006c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
13016c924f48SJed Brown     PetscFunctionReturn(0);
13026c924f48SJed Brown   }
130351f519a2SBarry Smith   for (i=0; i<n; i++) {
1304e32f2f54SBarry 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);
1305e32f2f54SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
130651f519a2SBarry Smith   }
1307b00a9115SJed Brown   ierr = PetscNew(&ilink);CHKERRQ(ierr);
1308a04f6461SBarry Smith   if (splitname) {
1309db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1310a04f6461SBarry Smith   } else {
1311785e854fSJed Brown     ierr = PetscMalloc1(3,&ilink->splitname);CHKERRQ(ierr);
1312a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
1313a04f6461SBarry Smith   }
13149df09d43SBarry 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 */
1315785e854fSJed Brown   ierr = PetscMalloc1(n,&ilink->fields);CHKERRQ(ierr);
13165a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
1317785e854fSJed Brown   ierr = PetscMalloc1(n,&ilink->fields_col);CHKERRQ(ierr);
13185d4c12cdSJungho Lee   ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr);
13192fa5cd67SKarl Rupp 
13205a9f2f41SSatish Balay   ilink->nfields = n;
13210298fd71SBarry Smith   ilink->next    = NULL;
1322ce94432eSBarry Smith   ierr           = KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);CHKERRQ(ierr);
1323422a814eSBarry Smith   ierr           = KSPSetErrorIfNotConverged(ilink->ksp,pc->erroriffailure);CHKERRQ(ierr);
132420252d06SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
13255a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
13269005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
132769a612a9SBarry Smith 
13288caf3d72SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);CHKERRQ(ierr);
13295a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
13300971522cSBarry Smith 
13310971522cSBarry Smith   if (!next) {
13325a9f2f41SSatish Balay     jac->head       = ilink;
13330298fd71SBarry Smith     ilink->previous = NULL;
13340971522cSBarry Smith   } else {
13350971522cSBarry Smith     while (next->next) {
13360971522cSBarry Smith       next = next->next;
13370971522cSBarry Smith     }
13385a9f2f41SSatish Balay     next->next      = ilink;
133951f519a2SBarry Smith     ilink->previous = next;
13400971522cSBarry Smith   }
13410971522cSBarry Smith   jac->nsplits++;
13420971522cSBarry Smith   PetscFunctionReturn(0);
13430971522cSBarry Smith }
13440971522cSBarry Smith 
1345e69d4d44SBarry Smith #undef __FUNCT__
1346e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
13471e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
1348e69d4d44SBarry Smith {
1349e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1350e69d4d44SBarry Smith   PetscErrorCode ierr;
1351e69d4d44SBarry Smith 
1352e69d4d44SBarry Smith   PetscFunctionBegin;
135334a3b412SBarry Smith   if (!jac->schur) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call KSPSetUp() or PCSetUp() before calling PCFieldSplitGetSubKSP()");
1354785e854fSJed Brown   ierr = PetscMalloc1(jac->nsplits,subksp);CHKERRQ(ierr);
1355e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
13562fa5cd67SKarl Rupp 
1357e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
135813e0d083SBarry Smith   if (n) *n = jac->nsplits;
1359e69d4d44SBarry Smith   PetscFunctionReturn(0);
1360e69d4d44SBarry Smith }
13610971522cSBarry Smith 
13620971522cSBarry Smith #undef __FUNCT__
136369a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
13641e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
13650971522cSBarry Smith {
13660971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
13670971522cSBarry Smith   PetscErrorCode    ierr;
13680971522cSBarry Smith   PetscInt          cnt   = 0;
13695a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
13700971522cSBarry Smith 
13710971522cSBarry Smith   PetscFunctionBegin;
1372785e854fSJed Brown   ierr = PetscMalloc1(jac->nsplits,subksp);CHKERRQ(ierr);
13735a9f2f41SSatish Balay   while (ilink) {
13745a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
13755a9f2f41SSatish Balay     ilink            = ilink->next;
13760971522cSBarry Smith   }
13775d480477SMatthew 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);
137813e0d083SBarry Smith   if (n) *n = jac->nsplits;
13790971522cSBarry Smith   PetscFunctionReturn(0);
13800971522cSBarry Smith }
13810971522cSBarry Smith 
1382704ba839SBarry Smith #undef __FUNCT__
13836dbb499eSCian Wilson #define __FUNCT__ "PCFieldSplitRestrictIS"
13846dbb499eSCian Wilson /*@C
13856dbb499eSCian Wilson     PCFieldSplitRestrictIS - Restricts the fieldsplit ISs to be within a given IS.
13866dbb499eSCian Wilson 
13876dbb499eSCian Wilson     Input Parameters:
13886dbb499eSCian Wilson +   pc  - the preconditioner context
13896dbb499eSCian Wilson +   is - the index set that defines the indices to which the fieldsplit is to be restricted
13906dbb499eSCian Wilson 
13916dbb499eSCian Wilson     Level: advanced
13926dbb499eSCian Wilson 
13936dbb499eSCian Wilson @*/
13946dbb499eSCian Wilson PetscErrorCode  PCFieldSplitRestrictIS(PC pc,IS isy)
13956dbb499eSCian Wilson {
13966dbb499eSCian Wilson   PetscErrorCode ierr;
13976dbb499eSCian Wilson 
13986dbb499eSCian Wilson   PetscFunctionBegin;
13996dbb499eSCian Wilson   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
14006dbb499eSCian Wilson   PetscValidHeaderSpecific(isy,IS_CLASSID,2);
14010246f55cSBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitRestrictIS_C",(PC,IS),(pc,isy));CHKERRQ(ierr);
14026dbb499eSCian Wilson   PetscFunctionReturn(0);
14036dbb499eSCian Wilson }
14046dbb499eSCian Wilson 
14056dbb499eSCian Wilson 
14066dbb499eSCian Wilson #undef __FUNCT__
14076dbb499eSCian Wilson #define __FUNCT__ "PCFieldSplitRestrictIS_FieldSplit"
14086dbb499eSCian Wilson static PetscErrorCode  PCFieldSplitRestrictIS_FieldSplit(PC pc, IS isy)
14096dbb499eSCian Wilson {
14106dbb499eSCian Wilson   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
14116dbb499eSCian Wilson   PetscErrorCode    ierr;
14126dbb499eSCian Wilson   PC_FieldSplitLink ilink = jac->head, next;
14136dbb499eSCian Wilson   PetscInt          localsize,size,sizez,i;
14146dbb499eSCian Wilson   const PetscInt    *ind, *indz;
14156dbb499eSCian Wilson   PetscInt          *indc, *indcz;
14166dbb499eSCian Wilson   PetscBool         flg;
14176dbb499eSCian Wilson 
14186dbb499eSCian Wilson   PetscFunctionBegin;
14196dbb499eSCian Wilson   ierr = ISGetLocalSize(isy,&localsize);CHKERRQ(ierr);
14206dbb499eSCian Wilson   ierr = MPI_Scan(&localsize,&size,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)isy));CHKERRQ(ierr);
14216dbb499eSCian Wilson   size -= localsize;
14226dbb499eSCian Wilson   while(ilink) {
14236dbb499eSCian Wilson     IS isrl,isr;
14241c7cfba8SBarry Smith     PC subpc;
14256dbb499eSCian Wilson     ierr          = ISEmbed(ilink->is, isy, PETSC_TRUE, &isrl);CHKERRQ(ierr);
14266dbb499eSCian Wilson     ierr          = ISGetLocalSize(isrl,&localsize);CHKERRQ(ierr);
14276dbb499eSCian Wilson     ierr          = PetscMalloc1(localsize,&indc);CHKERRQ(ierr);
14286dbb499eSCian Wilson     ierr          = ISGetIndices(isrl,&ind);CHKERRQ(ierr);
14296dbb499eSCian Wilson     ierr          = PetscMemcpy(indc,ind,localsize*sizeof(PetscInt));CHKERRQ(ierr);
14306dbb499eSCian Wilson     ierr          = ISRestoreIndices(isrl,&ind);CHKERRQ(ierr);
14316dbb499eSCian Wilson     ierr          = ISDestroy(&isrl);CHKERRQ(ierr);
14326dbb499eSCian Wilson     for (i=0; i<localsize; i++) *(indc+i) += size;
14336dbb499eSCian Wilson     ierr          = ISCreateGeneral(PetscObjectComm((PetscObject)isy),localsize,indc,PETSC_OWN_POINTER,&isr);CHKERRQ(ierr);
14346dbb499eSCian Wilson     ierr          = PetscObjectReference((PetscObject)isr);CHKERRQ(ierr);
14356dbb499eSCian Wilson     ierr          = ISDestroy(&ilink->is);CHKERRQ(ierr);
14366dbb499eSCian Wilson     ilink->is     = isr;
14376dbb499eSCian Wilson     ierr          = PetscObjectReference((PetscObject)isr);CHKERRQ(ierr);
14386dbb499eSCian Wilson     ierr          = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
14396dbb499eSCian Wilson     ilink->is_col = isr;
14406dbb499eSCian Wilson     ierr          = ISDestroy(&isr);CHKERRQ(ierr);
14416dbb499eSCian Wilson     ierr          = KSPGetPC(ilink->ksp, &subpc);CHKERRQ(ierr);
14426dbb499eSCian Wilson     ierr          = PetscObjectTypeCompare((PetscObject)subpc,PCFIELDSPLIT,&flg);CHKERRQ(ierr);
14436dbb499eSCian Wilson     if(flg) {
14446dbb499eSCian Wilson       IS iszl,isz;
14456dbb499eSCian Wilson       MPI_Comm comm;
14466dbb499eSCian Wilson       ierr   = ISGetLocalSize(ilink->is,&localsize);CHKERRQ(ierr);
14476dbb499eSCian Wilson       comm   = PetscObjectComm((PetscObject)ilink->is);
14486dbb499eSCian Wilson       ierr   = ISEmbed(isy, ilink->is, PETSC_TRUE, &iszl);CHKERRQ(ierr);
14496dbb499eSCian Wilson       ierr   = MPI_Scan(&localsize,&sizez,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr);
14506dbb499eSCian Wilson       sizez -= localsize;
14516dbb499eSCian Wilson       ierr   = ISGetLocalSize(iszl,&localsize);CHKERRQ(ierr);
14526dbb499eSCian Wilson       ierr   = PetscMalloc1(localsize,&indcz);CHKERRQ(ierr);
14536dbb499eSCian Wilson       ierr   = ISGetIndices(iszl,&indz);CHKERRQ(ierr);
14546dbb499eSCian Wilson       ierr   = PetscMemcpy(indcz,indz,localsize*sizeof(PetscInt));CHKERRQ(ierr);
14556dbb499eSCian Wilson       ierr   = ISRestoreIndices(iszl,&indz);CHKERRQ(ierr);
14566dbb499eSCian Wilson       ierr   = ISDestroy(&iszl);CHKERRQ(ierr);
14576dbb499eSCian Wilson       for (i=0; i<localsize; i++) *(indcz+i) += sizez;
14586dbb499eSCian Wilson       ierr   = ISCreateGeneral(comm,localsize,indcz,PETSC_OWN_POINTER,&isz);CHKERRQ(ierr);
14596dbb499eSCian Wilson       ierr   = PCFieldSplitRestrictIS(subpc,isz);CHKERRQ(ierr);
14606dbb499eSCian Wilson       ierr   = ISDestroy(&isz);CHKERRQ(ierr);
14616dbb499eSCian Wilson     }
14626dbb499eSCian Wilson     next = ilink->next;
14636dbb499eSCian Wilson     ilink = next;
14646dbb499eSCian Wilson   }
14656dbb499eSCian Wilson   jac->isrestrict = PETSC_TRUE;
14666dbb499eSCian Wilson   PetscFunctionReturn(0);
14676dbb499eSCian Wilson }
14686dbb499eSCian Wilson 
14696dbb499eSCian Wilson #undef __FUNCT__
1470704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
14711e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
1472704ba839SBarry Smith {
1473704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1474704ba839SBarry Smith   PetscErrorCode    ierr;
1475704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
1476704ba839SBarry Smith   char              prefix[128];
1477704ba839SBarry Smith 
1478704ba839SBarry Smith   PetscFunctionBegin;
14796c924f48SJed Brown   if (jac->splitdefined) {
14806c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
14816c924f48SJed Brown     PetscFunctionReturn(0);
14826c924f48SJed Brown   }
1483b00a9115SJed Brown   ierr = PetscNew(&ilink);CHKERRQ(ierr);
1484a04f6461SBarry Smith   if (splitname) {
1485db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1486a04f6461SBarry Smith   } else {
1487785e854fSJed Brown     ierr = PetscMalloc1(8,&ilink->splitname);CHKERRQ(ierr);
1488b5787286SJed Brown     ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr);
1489a04f6461SBarry Smith   }
14909df09d43SBarry 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 */
14911ab39975SBarry Smith   ierr          = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1492b5787286SJed Brown   ierr          = ISDestroy(&ilink->is);CHKERRQ(ierr);
1493b5787286SJed Brown   ilink->is     = is;
1494b5787286SJed Brown   ierr          = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1495b5787286SJed Brown   ierr          = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
1496b5787286SJed Brown   ilink->is_col = is;
14970298fd71SBarry Smith   ilink->next   = NULL;
1498ce94432eSBarry Smith   ierr          = KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);CHKERRQ(ierr);
1499422a814eSBarry Smith   ierr          = KSPSetErrorIfNotConverged(ilink->ksp,pc->erroriffailure);CHKERRQ(ierr);
150020252d06SBarry Smith   ierr          = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1501704ba839SBarry Smith   ierr          = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
15029005cf84SBarry Smith   ierr          = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
1503704ba839SBarry Smith 
15048caf3d72SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);CHKERRQ(ierr);
1505704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
1506704ba839SBarry Smith 
1507704ba839SBarry Smith   if (!next) {
1508704ba839SBarry Smith     jac->head       = ilink;
15090298fd71SBarry Smith     ilink->previous = NULL;
1510704ba839SBarry Smith   } else {
1511704ba839SBarry Smith     while (next->next) {
1512704ba839SBarry Smith       next = next->next;
1513704ba839SBarry Smith     }
1514704ba839SBarry Smith     next->next      = ilink;
1515704ba839SBarry Smith     ilink->previous = next;
1516704ba839SBarry Smith   }
1517704ba839SBarry Smith   jac->nsplits++;
1518704ba839SBarry Smith   PetscFunctionReturn(0);
1519704ba839SBarry Smith }
1520704ba839SBarry Smith 
15210971522cSBarry Smith #undef __FUNCT__
15220971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
15230971522cSBarry Smith /*@
15240971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
15250971522cSBarry Smith 
1526ad4df100SBarry Smith     Logically Collective on PC
15270971522cSBarry Smith 
15280971522cSBarry Smith     Input Parameters:
15290971522cSBarry Smith +   pc  - the preconditioner context
15300298fd71SBarry Smith .   splitname - name of this split, if NULL the number of the split is used
15310971522cSBarry Smith .   n - the number of fields in this split
1532db4c96c1SJed Brown -   fields - the fields in this split
15330971522cSBarry Smith 
15340971522cSBarry Smith     Level: intermediate
15350971522cSBarry Smith 
1536d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
1537d32f9abdSBarry Smith 
15387287d2a3SDmitry Karpeev      The PCFieldSplitSetFields() is for defining fields as strided blocks. For example, if the block
1539d32f9abdSBarry 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
1540d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1541d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
1542d32f9abdSBarry Smith 
1543db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
1544db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1545db4c96c1SJed Brown 
15465d4c12cdSJungho Lee      Developer Note: This routine does not actually create the IS representing the split, that is delayed
15475d4c12cdSJungho Lee      until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be
15485d4c12cdSJungho Lee      available when this routine is called.
15495d4c12cdSJungho Lee 
1550d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
15510971522cSBarry Smith 
15520971522cSBarry Smith @*/
15535d4c12cdSJungho Lee PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
15540971522cSBarry Smith {
15554ac538c5SBarry Smith   PetscErrorCode ierr;
15560971522cSBarry Smith 
15570971522cSBarry Smith   PetscFunctionBegin;
15580700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1559db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1560ce94432eSBarry Smith   if (n < 1) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1561db4c96c1SJed Brown   PetscValidIntPointer(fields,3);
15620246f55cSBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt*,const PetscInt*),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr);
15630971522cSBarry Smith   PetscFunctionReturn(0);
15640971522cSBarry Smith }
15650971522cSBarry Smith 
15660971522cSBarry Smith #undef __FUNCT__
1567c84da90fSDmitry Karpeev #define __FUNCT__ "PCFieldSplitSetDiagUseAmat"
1568c84da90fSDmitry Karpeev /*@
1569c84da90fSDmitry Karpeev     PCFieldSplitSetDiagUseAmat - set flag indicating whether to extract diagonal blocks from Amat (rather than Pmat)
1570c84da90fSDmitry Karpeev 
1571c84da90fSDmitry Karpeev     Logically Collective on PC
1572c84da90fSDmitry Karpeev 
1573c84da90fSDmitry Karpeev     Input Parameters:
1574c84da90fSDmitry Karpeev +   pc  - the preconditioner object
1575c84da90fSDmitry Karpeev -   flg - boolean flag indicating whether or not to use Amat to extract the diagonal blocks from
1576c84da90fSDmitry Karpeev 
15779506b623SDmitry Karpeev     Options Database:
15789506b623SDmitry Karpeev .     -pc_fieldsplit_diag_use_amat
1579c84da90fSDmitry Karpeev 
1580c84da90fSDmitry Karpeev     Level: intermediate
1581c84da90fSDmitry Karpeev 
1582c84da90fSDmitry Karpeev .seealso: PCFieldSplitGetDiagUseAmat(), PCFieldSplitSetOffDiagUseAmat(), PCFIELDSPLIT
1583c84da90fSDmitry Karpeev 
1584c84da90fSDmitry Karpeev @*/
1585c84da90fSDmitry Karpeev PetscErrorCode  PCFieldSplitSetDiagUseAmat(PC pc,PetscBool flg)
1586c84da90fSDmitry Karpeev {
1587c84da90fSDmitry Karpeev   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1588c84da90fSDmitry Karpeev   PetscBool      isfs;
1589c84da90fSDmitry Karpeev   PetscErrorCode ierr;
1590c84da90fSDmitry Karpeev 
1591c84da90fSDmitry Karpeev   PetscFunctionBegin;
1592c84da90fSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1593c84da90fSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
1594c84da90fSDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1595c84da90fSDmitry Karpeev   jac->diag_use_amat = flg;
1596c84da90fSDmitry Karpeev   PetscFunctionReturn(0);
1597c84da90fSDmitry Karpeev }
1598c84da90fSDmitry Karpeev 
1599c84da90fSDmitry Karpeev #undef __FUNCT__
1600c84da90fSDmitry Karpeev #define __FUNCT__ "PCFieldSplitGetDiagUseAmat"
1601c84da90fSDmitry Karpeev /*@
1602c84da90fSDmitry Karpeev     PCFieldSplitGetDiagUseAmat - get the flag indicating whether to extract diagonal blocks from Amat (rather than Pmat)
1603c84da90fSDmitry Karpeev 
1604c84da90fSDmitry Karpeev     Logically Collective on PC
1605c84da90fSDmitry Karpeev 
1606c84da90fSDmitry Karpeev     Input Parameters:
1607c84da90fSDmitry Karpeev .   pc  - the preconditioner object
1608c84da90fSDmitry Karpeev 
1609c84da90fSDmitry Karpeev     Output Parameters:
1610c84da90fSDmitry Karpeev .   flg - boolean flag indicating whether or not to use Amat to extract the diagonal blocks from
1611c84da90fSDmitry Karpeev 
1612c84da90fSDmitry Karpeev 
1613c84da90fSDmitry Karpeev     Level: intermediate
1614c84da90fSDmitry Karpeev 
1615c84da90fSDmitry Karpeev .seealso: PCFieldSplitSetDiagUseAmat(), PCFieldSplitGetOffDiagUseAmat(), PCFIELDSPLIT
1616c84da90fSDmitry Karpeev 
1617c84da90fSDmitry Karpeev @*/
1618c84da90fSDmitry Karpeev PetscErrorCode  PCFieldSplitGetDiagUseAmat(PC pc,PetscBool *flg)
1619c84da90fSDmitry Karpeev {
1620c84da90fSDmitry Karpeev   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1621c84da90fSDmitry Karpeev   PetscBool      isfs;
1622c84da90fSDmitry Karpeev   PetscErrorCode ierr;
1623c84da90fSDmitry Karpeev 
1624c84da90fSDmitry Karpeev   PetscFunctionBegin;
1625c84da90fSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1626c84da90fSDmitry Karpeev   PetscValidPointer(flg,2);
1627c84da90fSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
1628c84da90fSDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1629c84da90fSDmitry Karpeev   *flg = jac->diag_use_amat;
1630c84da90fSDmitry Karpeev   PetscFunctionReturn(0);
1631c84da90fSDmitry Karpeev }
1632c84da90fSDmitry Karpeev 
1633c84da90fSDmitry Karpeev #undef __FUNCT__
1634c84da90fSDmitry Karpeev #define __FUNCT__ "PCFieldSplitSetOffDiagUseAmat"
1635c84da90fSDmitry Karpeev /*@
1636c84da90fSDmitry Karpeev     PCFieldSplitSetOffDiagUseAmat - set flag indicating whether to extract off-diagonal blocks from Amat (rather than Pmat)
1637c84da90fSDmitry Karpeev 
1638c84da90fSDmitry Karpeev     Logically Collective on PC
1639c84da90fSDmitry Karpeev 
1640c84da90fSDmitry Karpeev     Input Parameters:
1641c84da90fSDmitry Karpeev +   pc  - the preconditioner object
1642c84da90fSDmitry Karpeev -   flg - boolean flag indicating whether or not to use Amat to extract the off-diagonal blocks from
1643c84da90fSDmitry Karpeev 
16449506b623SDmitry Karpeev     Options Database:
16459506b623SDmitry Karpeev .     -pc_fieldsplit_off_diag_use_amat
1646c84da90fSDmitry Karpeev 
1647c84da90fSDmitry Karpeev     Level: intermediate
1648c84da90fSDmitry Karpeev 
1649c84da90fSDmitry Karpeev .seealso: PCFieldSplitGetOffDiagUseAmat(), PCFieldSplitSetDiagUseAmat(), PCFIELDSPLIT
1650c84da90fSDmitry Karpeev 
1651c84da90fSDmitry Karpeev @*/
1652c84da90fSDmitry Karpeev PetscErrorCode  PCFieldSplitSetOffDiagUseAmat(PC pc,PetscBool flg)
1653c84da90fSDmitry Karpeev {
1654c84da90fSDmitry Karpeev   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1655c84da90fSDmitry Karpeev   PetscBool      isfs;
1656c84da90fSDmitry Karpeev   PetscErrorCode ierr;
1657c84da90fSDmitry Karpeev 
1658c84da90fSDmitry Karpeev   PetscFunctionBegin;
1659c84da90fSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1660c84da90fSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
1661c84da90fSDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1662c84da90fSDmitry Karpeev   jac->offdiag_use_amat = flg;
1663c84da90fSDmitry Karpeev   PetscFunctionReturn(0);
1664c84da90fSDmitry Karpeev }
1665c84da90fSDmitry Karpeev 
1666c84da90fSDmitry Karpeev #undef __FUNCT__
1667c84da90fSDmitry Karpeev #define __FUNCT__ "PCFieldSplitGetOffDiagUseAmat"
1668c84da90fSDmitry Karpeev /*@
1669c84da90fSDmitry Karpeev     PCFieldSplitGetOffDiagUseAmat - get the flag indicating whether to extract off-diagonal blocks from Amat (rather than Pmat)
1670c84da90fSDmitry Karpeev 
1671c84da90fSDmitry Karpeev     Logically Collective on PC
1672c84da90fSDmitry Karpeev 
1673c84da90fSDmitry Karpeev     Input Parameters:
1674c84da90fSDmitry Karpeev .   pc  - the preconditioner object
1675c84da90fSDmitry Karpeev 
1676c84da90fSDmitry Karpeev     Output Parameters:
1677c84da90fSDmitry Karpeev .   flg - boolean flag indicating whether or not to use Amat to extract the off-diagonal blocks from
1678c84da90fSDmitry Karpeev 
1679c84da90fSDmitry Karpeev 
1680c84da90fSDmitry Karpeev     Level: intermediate
1681c84da90fSDmitry Karpeev 
1682c84da90fSDmitry Karpeev .seealso: PCFieldSplitSetOffDiagUseAmat(), PCFieldSplitGetDiagUseAmat(), PCFIELDSPLIT
1683c84da90fSDmitry Karpeev 
1684c84da90fSDmitry Karpeev @*/
1685c84da90fSDmitry Karpeev PetscErrorCode  PCFieldSplitGetOffDiagUseAmat(PC pc,PetscBool *flg)
1686c84da90fSDmitry Karpeev {
1687c84da90fSDmitry Karpeev   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1688c84da90fSDmitry Karpeev   PetscBool      isfs;
1689c84da90fSDmitry Karpeev   PetscErrorCode ierr;
1690c84da90fSDmitry Karpeev 
1691c84da90fSDmitry Karpeev   PetscFunctionBegin;
1692c84da90fSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1693c84da90fSDmitry Karpeev   PetscValidPointer(flg,2);
1694c84da90fSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
1695c84da90fSDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"PC not of type %s",PCFIELDSPLIT);
1696c84da90fSDmitry Karpeev   *flg = jac->offdiag_use_amat;
1697c84da90fSDmitry Karpeev   PetscFunctionReturn(0);
1698c84da90fSDmitry Karpeev }
1699c84da90fSDmitry Karpeev 
1700c84da90fSDmitry Karpeev 
1701c84da90fSDmitry Karpeev 
1702c84da90fSDmitry Karpeev #undef __FUNCT__
1703704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
1704218d4943SBarry Smith /*@C
1705704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
1706704ba839SBarry Smith 
1707ad4df100SBarry Smith     Logically Collective on PC
1708704ba839SBarry Smith 
1709704ba839SBarry Smith     Input Parameters:
1710704ba839SBarry Smith +   pc  - the preconditioner context
17110298fd71SBarry Smith .   splitname - name of this split, if NULL the number of the split is used
1712db4c96c1SJed Brown -   is - the index set that defines the vector elements in this field
1713704ba839SBarry Smith 
1714d32f9abdSBarry Smith 
1715a6ffb8dbSJed Brown     Notes:
1716a6ffb8dbSJed Brown     Use PCFieldSplitSetFields(), for fields defined by strided types.
1717a6ffb8dbSJed Brown 
1718db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
1719db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1720d32f9abdSBarry Smith 
1721704ba839SBarry Smith     Level: intermediate
1722704ba839SBarry Smith 
1723704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
1724704ba839SBarry Smith 
1725704ba839SBarry Smith @*/
17267087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1727704ba839SBarry Smith {
17284ac538c5SBarry Smith   PetscErrorCode ierr;
1729704ba839SBarry Smith 
1730704ba839SBarry Smith   PetscFunctionBegin;
17310700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
17327b62db95SJungho Lee   if (splitname) PetscValidCharPointer(splitname,2);
1733db4c96c1SJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,3);
1734ccc738f7SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr);
1735704ba839SBarry Smith   PetscFunctionReturn(0);
1736704ba839SBarry Smith }
1737704ba839SBarry Smith 
1738704ba839SBarry Smith #undef __FUNCT__
173957a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS"
174057a9adfeSMatthew G Knepley /*@
174157a9adfeSMatthew G Knepley     PCFieldSplitGetIS - Retrieves the elements for a field as an IS
174257a9adfeSMatthew G Knepley 
174357a9adfeSMatthew G Knepley     Logically Collective on PC
174457a9adfeSMatthew G Knepley 
174557a9adfeSMatthew G Knepley     Input Parameters:
174657a9adfeSMatthew G Knepley +   pc  - the preconditioner context
174757a9adfeSMatthew G Knepley -   splitname - name of this split
174857a9adfeSMatthew G Knepley 
174957a9adfeSMatthew G Knepley     Output Parameter:
17500298fd71SBarry Smith -   is - the index set that defines the vector elements in this field, or NULL if the field is not found
175157a9adfeSMatthew G Knepley 
175257a9adfeSMatthew G Knepley     Level: intermediate
175357a9adfeSMatthew G Knepley 
175457a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS()
175557a9adfeSMatthew G Knepley 
175657a9adfeSMatthew G Knepley @*/
175757a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
175857a9adfeSMatthew G Knepley {
175957a9adfeSMatthew G Knepley   PetscErrorCode ierr;
176057a9adfeSMatthew G Knepley 
176157a9adfeSMatthew G Knepley   PetscFunctionBegin;
176257a9adfeSMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
176357a9adfeSMatthew G Knepley   PetscValidCharPointer(splitname,2);
176457a9adfeSMatthew G Knepley   PetscValidPointer(is,3);
176557a9adfeSMatthew G Knepley   {
176657a9adfeSMatthew G Knepley     PC_FieldSplit     *jac  = (PC_FieldSplit*) pc->data;
176757a9adfeSMatthew G Knepley     PC_FieldSplitLink ilink = jac->head;
176857a9adfeSMatthew G Knepley     PetscBool         found;
176957a9adfeSMatthew G Knepley 
17700298fd71SBarry Smith     *is = NULL;
177157a9adfeSMatthew G Knepley     while (ilink) {
177257a9adfeSMatthew G Knepley       ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr);
177357a9adfeSMatthew G Knepley       if (found) {
177457a9adfeSMatthew G Knepley         *is = ilink->is;
177557a9adfeSMatthew G Knepley         break;
177657a9adfeSMatthew G Knepley       }
177757a9adfeSMatthew G Knepley       ilink = ilink->next;
177857a9adfeSMatthew G Knepley     }
177957a9adfeSMatthew G Knepley   }
178057a9adfeSMatthew G Knepley   PetscFunctionReturn(0);
178157a9adfeSMatthew G Knepley }
178257a9adfeSMatthew G Knepley 
178357a9adfeSMatthew G Knepley #undef __FUNCT__
178451f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
178551f519a2SBarry Smith /*@
178651f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
178751f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
178851f519a2SBarry Smith 
1789ad4df100SBarry Smith     Logically Collective on PC
179051f519a2SBarry Smith 
179151f519a2SBarry Smith     Input Parameters:
179251f519a2SBarry Smith +   pc  - the preconditioner context
179351f519a2SBarry Smith -   bs - the block size
179451f519a2SBarry Smith 
179551f519a2SBarry Smith     Level: intermediate
179651f519a2SBarry Smith 
179751f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
179851f519a2SBarry Smith 
179951f519a2SBarry Smith @*/
18007087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
180151f519a2SBarry Smith {
18024ac538c5SBarry Smith   PetscErrorCode ierr;
180351f519a2SBarry Smith 
180451f519a2SBarry Smith   PetscFunctionBegin;
18050700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1806c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,bs,2);
18074ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr);
180851f519a2SBarry Smith   PetscFunctionReturn(0);
180951f519a2SBarry Smith }
181051f519a2SBarry Smith 
181151f519a2SBarry Smith #undef __FUNCT__
181269a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
18130971522cSBarry Smith /*@C
181469a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
18150971522cSBarry Smith 
181669a612a9SBarry Smith    Collective on KSP
18170971522cSBarry Smith 
18180971522cSBarry Smith    Input Parameter:
18190971522cSBarry Smith .  pc - the preconditioner context
18200971522cSBarry Smith 
18210971522cSBarry Smith    Output Parameters:
182213e0d083SBarry Smith +  n - the number of splits
18233111de3cSBarry Smith -  subksp - the array of KSP contexts
18240971522cSBarry Smith 
18250971522cSBarry Smith    Note:
18269a4f7e47SBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs is to be freed by the user with PetscFree()
1827d32f9abdSBarry Smith    (not the KSP just the array that contains them).
18280971522cSBarry Smith 
182969a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
18300971522cSBarry Smith 
1831196cc216SBarry Smith    Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs.
1832*2bf68e3eSBarry Smith       You can call PCFieldSplitGetSubKSP(pc,n,PETSC_NULL_KSP,ierr) to determine how large the
1833196cc216SBarry Smith       KSP array must be.
1834196cc216SBarry Smith 
1835196cc216SBarry Smith 
18360971522cSBarry Smith    Level: advanced
18370971522cSBarry Smith 
18380971522cSBarry Smith .seealso: PCFIELDSPLIT
18390971522cSBarry Smith @*/
18407087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
18410971522cSBarry Smith {
18424ac538c5SBarry Smith   PetscErrorCode ierr;
18430971522cSBarry Smith 
18440971522cSBarry Smith   PetscFunctionBegin;
18450700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
184613e0d083SBarry Smith   if (n) PetscValidIntPointer(n,2);
18474ac538c5SBarry Smith   ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
18480971522cSBarry Smith   PetscFunctionReturn(0);
18490971522cSBarry Smith }
18500971522cSBarry Smith 
1851e69d4d44SBarry Smith #undef __FUNCT__
185229f8a81cSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurPre"
1853e69d4d44SBarry Smith /*@
185453f2277eSBarry Smith     PCFieldSplitSetSchurPre -  Indicates what operator is used to construct the preconditioner for the Schur complement.
1855a04f6461SBarry Smith       A11 matrix. Otherwise no preconditioner is used.
1856e69d4d44SBarry Smith 
1857e69d4d44SBarry Smith     Collective on PC
1858e69d4d44SBarry Smith 
1859e69d4d44SBarry Smith     Input Parameters:
1860e69d4d44SBarry Smith +   pc      - the preconditioner context
18616fdd48a9SBarry 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
18626fdd48a9SBarry Smith               PC_FIELDSPLIT_SCHUR_PRE_SELFP, and PC_FIELDSPLIT_SCHUR_PRE_FULL
18630298fd71SBarry Smith -   userpre - matrix to use for preconditioning, or NULL
1864084e4875SJed Brown 
1865e69d4d44SBarry Smith     Options Database:
18666fdd48a9SBarry Smith .     -pc_fieldsplit_schur_precondition <self,selfp,user,a11,full> - default is a11. See notes for meaning of various arguments
1867e69d4d44SBarry Smith 
1868fd1303e9SJungho Lee     Notes:
1869fd1303e9SJungho Lee $    If ptype is
1870a6a584a2SBarry Smith $        a11 then the preconditioner for the Schur complement is generated from the block diagonal part of the preconditioner
187153f2277eSBarry Smith $             matrix associated with the Schur complement (i.e. A11), not the Schur complement matrix
1872a6a584a2SBarry Smith $        self the preconditioner for the Schur complement is generated from the symbolic representation of the Schur complement matrix:
1873a6a584a2SBarry Smith $             The only preconditioner that currently works with this symbolic respresentation matrix object is the PCLSC
1874fd1303e9SJungho Lee $             preconditioner
18756fdd48a9SBarry Smith $        user then the preconditioner for the Schur complement is generated from the user provided matrix (pre argument
18766fdd48a9SBarry Smith $             to this function).
1877a6a584a2SBarry Smith $        selfp then the preconditioning for the Schur complement is generated from an explicitly-assembled approximation Sp = A11 - A10 inv(diag(A00)) A01
18781d71051fSDmitry Karpeev $             This is only a good preconditioner when diag(A00) is a good preconditioner for A00. Optionally, A00 can be
18796fdd48a9SBarry Smith $             lumped before extracting the diagonal using the additional option -fieldsplit_1_mat_schur_complement_ainv_type lump
18806fdd48a9SBarry 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)
18816fdd48a9SBarry Smith $             useful mostly as a test that the Schur complement approach can work for your problem
1882fd1303e9SJungho Lee 
1883e87fab1bSBarry Smith      When solving a saddle point problem, where the A11 block is identically zero, using a11 as the ptype only makes sense
1884fd1303e9SJungho Lee     with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
1885a7476a74SDmitry Karpeev     -fieldsplit_1_pc_type lsc which uses the least squares commutator to compute a preconditioner for the Schur complement.
1886fd1303e9SJungho Lee 
1887e69d4d44SBarry Smith     Level: intermediate
1888e69d4d44SBarry Smith 
18891d71051fSDmitry Karpeev .seealso: PCFieldSplitGetSchurPre(), PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType,
18901d71051fSDmitry Karpeev           MatSchurComplementSetAinvType(), PCLSC
1891e69d4d44SBarry Smith 
1892e69d4d44SBarry Smith @*/
189329f8a81cSJed Brown PetscErrorCode PCFieldSplitSetSchurPre(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1894e69d4d44SBarry Smith {
18954ac538c5SBarry Smith   PetscErrorCode ierr;
1896e69d4d44SBarry Smith 
1897e69d4d44SBarry Smith   PetscFunctionBegin;
18980700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
189929f8a81cSJed Brown   ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurPre_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr);
1900e69d4d44SBarry Smith   PetscFunctionReturn(0);
1901e69d4d44SBarry Smith }
190229f8a81cSJed Brown PetscErrorCode PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) {return PCFieldSplitSetSchurPre(pc,ptype,pre);} /* Deprecated name */
1903e69d4d44SBarry Smith 
1904e69d4d44SBarry Smith #undef __FUNCT__
190537a82bf0SJed Brown #define __FUNCT__ "PCFieldSplitGetSchurPre"
190637a82bf0SJed Brown /*@
190737a82bf0SJed Brown     PCFieldSplitGetSchurPre - For Schur complement fieldsplit, determine how the Schur complement will be
190829f8a81cSJed Brown     preconditioned.  See PCFieldSplitSetSchurPre() for details.
190937a82bf0SJed Brown 
191037a82bf0SJed Brown     Logically Collective on PC
191137a82bf0SJed Brown 
191237a82bf0SJed Brown     Input Parameters:
191337a82bf0SJed Brown .   pc      - the preconditioner context
191437a82bf0SJed Brown 
191537a82bf0SJed Brown     Output Parameters:
191637a82bf0SJed 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
191737a82bf0SJed Brown -   userpre - matrix to use for preconditioning (with PC_FIELDSPLIT_PRE_USER), or NULL
191837a82bf0SJed Brown 
191937a82bf0SJed Brown     Level: intermediate
192037a82bf0SJed Brown 
192129f8a81cSJed Brown .seealso: PCFieldSplitSetSchurPre(), PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC
192237a82bf0SJed Brown 
192337a82bf0SJed Brown @*/
192437a82bf0SJed Brown PetscErrorCode PCFieldSplitGetSchurPre(PC pc,PCFieldSplitSchurPreType *ptype,Mat *pre)
192537a82bf0SJed Brown {
192637a82bf0SJed Brown   PetscErrorCode ierr;
192737a82bf0SJed Brown 
192837a82bf0SJed Brown   PetscFunctionBegin;
192937a82bf0SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
193037a82bf0SJed Brown   ierr = PetscUseMethod(pc,"PCFieldSplitGetSchurPre_C",(PC,PCFieldSplitSchurPreType*,Mat*),(pc,ptype,pre));CHKERRQ(ierr);
1931e69d4d44SBarry Smith   PetscFunctionReturn(0);
1932e69d4d44SBarry Smith }
1933e69d4d44SBarry Smith 
1934e69d4d44SBarry Smith #undef __FUNCT__
1935470b340bSDmitry Karpeev #define __FUNCT__ "PCFieldSplitSchurGetS"
193645e7fc46SDmitry Karpeev /*@
1937470b340bSDmitry Karpeev     PCFieldSplitSchurGetS -  extract the MatSchurComplement object used by this PC in case it needs to be configured separately
193845e7fc46SDmitry Karpeev 
193945e7fc46SDmitry Karpeev     Not collective
194045e7fc46SDmitry Karpeev 
194145e7fc46SDmitry Karpeev     Input Parameter:
194245e7fc46SDmitry Karpeev .   pc      - the preconditioner context
194345e7fc46SDmitry Karpeev 
194445e7fc46SDmitry Karpeev     Output Parameter:
1945470b340bSDmitry Karpeev .   S       - the Schur complement matrix
194645e7fc46SDmitry Karpeev 
1947470b340bSDmitry Karpeev     Notes:
1948470b340bSDmitry Karpeev     This matrix should not be destroyed using MatDestroy(); rather, use PCFieldSplitSchurRestoreS().
194945e7fc46SDmitry Karpeev 
195045e7fc46SDmitry Karpeev     Level: advanced
195145e7fc46SDmitry Karpeev 
195229f8a81cSJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSchurPreType, PCFieldSplitSetSchurPre(), MatSchurComplement, PCFieldSplitSchurRestoreS()
195345e7fc46SDmitry Karpeev 
195445e7fc46SDmitry Karpeev @*/
1955470b340bSDmitry Karpeev PetscErrorCode  PCFieldSplitSchurGetS(PC pc,Mat *S)
195645e7fc46SDmitry Karpeev {
195745e7fc46SDmitry Karpeev   PetscErrorCode ierr;
195845e7fc46SDmitry Karpeev   const char*    t;
195945e7fc46SDmitry Karpeev   PetscBool      isfs;
196045e7fc46SDmitry Karpeev   PC_FieldSplit  *jac;
196145e7fc46SDmitry Karpeev 
196245e7fc46SDmitry Karpeev   PetscFunctionBegin;
196345e7fc46SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
196445e7fc46SDmitry Karpeev   ierr = PetscObjectGetType((PetscObject)pc,&t);CHKERRQ(ierr);
196545e7fc46SDmitry Karpeev   ierr = PetscStrcmp(t,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
196645e7fc46SDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PC of type PCFIELDSPLIT, got %s instead",t);
196745e7fc46SDmitry Karpeev   jac = (PC_FieldSplit*)pc->data;
1968470b340bSDmitry 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);
1969470b340bSDmitry Karpeev   if (S) *S = jac->schur;
197045e7fc46SDmitry Karpeev   PetscFunctionReturn(0);
197145e7fc46SDmitry Karpeev }
197245e7fc46SDmitry Karpeev 
197345e7fc46SDmitry Karpeev #undef __FUNCT__
1974470b340bSDmitry Karpeev #define __FUNCT__ "PCFieldSplitSchurRestoreS"
1975470b340bSDmitry Karpeev /*@
1976470b340bSDmitry Karpeev     PCFieldSplitSchurRestoreS -  restores the MatSchurComplement object used by this PC
1977470b340bSDmitry Karpeev 
1978470b340bSDmitry Karpeev     Not collective
1979470b340bSDmitry Karpeev 
1980470b340bSDmitry Karpeev     Input Parameters:
1981470b340bSDmitry Karpeev +   pc      - the preconditioner context
1982470b340bSDmitry Karpeev .   S       - the Schur complement matrix
1983470b340bSDmitry Karpeev 
1984470b340bSDmitry Karpeev     Level: advanced
1985470b340bSDmitry Karpeev 
198629f8a81cSJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSchurPreType, PCFieldSplitSetSchurPre(), MatSchurComplement, PCFieldSplitSchurGetS()
1987470b340bSDmitry Karpeev 
1988470b340bSDmitry Karpeev @*/
1989bca69d2bSDmitry Karpeev PetscErrorCode  PCFieldSplitSchurRestoreS(PC pc,Mat *S)
1990470b340bSDmitry Karpeev {
1991470b340bSDmitry Karpeev   PetscErrorCode ierr;
1992470b340bSDmitry Karpeev   const char*    t;
1993470b340bSDmitry Karpeev   PetscBool      isfs;
1994470b340bSDmitry Karpeev   PC_FieldSplit  *jac;
1995470b340bSDmitry Karpeev 
1996470b340bSDmitry Karpeev   PetscFunctionBegin;
1997470b340bSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1998470b340bSDmitry Karpeev   ierr = PetscObjectGetType((PetscObject)pc,&t);CHKERRQ(ierr);
1999470b340bSDmitry Karpeev   ierr = PetscStrcmp(t,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
2000470b340bSDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PC of type PCFIELDSPLIT, got %s instead",t);
2001470b340bSDmitry Karpeev   jac = (PC_FieldSplit*)pc->data;
2002470b340bSDmitry 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);
2003bca69d2bSDmitry Karpeev   if (!S || *S != jac->schur) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"MatSchurComplement restored is not the same as gotten");
2004470b340bSDmitry Karpeev   PetscFunctionReturn(0);
2005470b340bSDmitry Karpeev }
2006470b340bSDmitry Karpeev 
2007470b340bSDmitry Karpeev 
2008470b340bSDmitry Karpeev #undef __FUNCT__
200929f8a81cSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurPre_FieldSplit"
201029f8a81cSJed Brown static PetscErrorCode  PCFieldSplitSetSchurPre_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
2011e69d4d44SBarry Smith {
2012e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
2013084e4875SJed Brown   PetscErrorCode ierr;
2014e69d4d44SBarry Smith 
2015e69d4d44SBarry Smith   PetscFunctionBegin;
2016084e4875SJed Brown   jac->schurpre = ptype;
2017a7476a74SDmitry Karpeev   if (ptype == PC_FIELDSPLIT_SCHUR_PRE_USER && pre) {
20186bf464f9SBarry Smith     ierr            = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
2019084e4875SJed Brown     jac->schur_user = pre;
2020084e4875SJed Brown     ierr            = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
2021084e4875SJed Brown   }
2022e69d4d44SBarry Smith   PetscFunctionReturn(0);
2023e69d4d44SBarry Smith }
2024e69d4d44SBarry Smith 
202530ad9308SMatthew Knepley #undef __FUNCT__
202637a82bf0SJed Brown #define __FUNCT__ "PCFieldSplitGetSchurPre_FieldSplit"
202737a82bf0SJed Brown static PetscErrorCode  PCFieldSplitGetSchurPre_FieldSplit(PC pc,PCFieldSplitSchurPreType *ptype,Mat *pre)
202837a82bf0SJed Brown {
202937a82bf0SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
203037a82bf0SJed Brown 
203137a82bf0SJed Brown   PetscFunctionBegin;
203237a82bf0SJed Brown   *ptype = jac->schurpre;
203337a82bf0SJed Brown   *pre   = jac->schur_user;
203437a82bf0SJed Brown   PetscFunctionReturn(0);
203537a82bf0SJed Brown }
203637a82bf0SJed Brown 
203737a82bf0SJed Brown #undef __FUNCT__
2038c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType"
2039ab1df9f5SJed Brown /*@
2040c9c6ffaaSJed Brown     PCFieldSplitSetSchurFactType -  sets which blocks of the approximate block factorization to retain
2041ab1df9f5SJed Brown 
2042ab1df9f5SJed Brown     Collective on PC
2043ab1df9f5SJed Brown 
2044ab1df9f5SJed Brown     Input Parameters:
2045ab1df9f5SJed Brown +   pc  - the preconditioner context
2046c9c6ffaaSJed Brown -   ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default
2047ab1df9f5SJed Brown 
2048ab1df9f5SJed Brown     Options Database:
2049c9c6ffaaSJed Brown .     -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full
2050ab1df9f5SJed Brown 
2051ab1df9f5SJed Brown 
2052ab1df9f5SJed Brown     Level: intermediate
2053ab1df9f5SJed Brown 
2054ab1df9f5SJed Brown     Notes:
2055ab1df9f5SJed Brown     The FULL factorization is
2056ab1df9f5SJed Brown 
2057ab1df9f5SJed Brown $   (A   B)  = (1       0) (A   0) (1  Ainv*B)
2058ab1df9f5SJed Brown $   (C   D)    (C*Ainv  1) (0   S) (0     1  )
2059ab1df9f5SJed Brown 
20606be4592eSBarry 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,
20616be4592eSBarry 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).
2062ab1df9f5SJed Brown 
20636be4592eSBarry 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
20646be4592eSBarry 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
20656be4592eSBarry 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,
20666be4592eSBarry 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.
2067ab1df9f5SJed Brown 
20686be4592eSBarry 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
20696be4592eSBarry 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).
2070ab1df9f5SJed Brown 
2071ab1df9f5SJed Brown     References:
207296a0c994SBarry Smith +   1. - Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000).
207396a0c994SBarry Smith -   2. - Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001).
2074ab1df9f5SJed Brown 
2075ab1df9f5SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType
2076ab1df9f5SJed Brown @*/
2077c9c6ffaaSJed Brown PetscErrorCode  PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype)
2078ab1df9f5SJed Brown {
2079ab1df9f5SJed Brown   PetscErrorCode ierr;
2080ab1df9f5SJed Brown 
2081ab1df9f5SJed Brown   PetscFunctionBegin;
2082ab1df9f5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2083c9c6ffaaSJed Brown   ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr);
2084ab1df9f5SJed Brown   PetscFunctionReturn(0);
2085ab1df9f5SJed Brown }
2086ab1df9f5SJed Brown 
2087ab1df9f5SJed Brown #undef __FUNCT__
2088c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType_FieldSplit"
20891e6b0712SBarry Smith static PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype)
2090ab1df9f5SJed Brown {
2091ab1df9f5SJed Brown   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
2092ab1df9f5SJed Brown 
2093ab1df9f5SJed Brown   PetscFunctionBegin;
2094ab1df9f5SJed Brown   jac->schurfactorization = ftype;
2095ab1df9f5SJed Brown   PetscFunctionReturn(0);
2096ab1df9f5SJed Brown }
2097ab1df9f5SJed Brown 
2098ab1df9f5SJed Brown #undef __FUNCT__
209930ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
210030ad9308SMatthew Knepley /*@C
21018c03b21aSDmitry Karpeev    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
210230ad9308SMatthew Knepley 
210330ad9308SMatthew Knepley    Collective on KSP
210430ad9308SMatthew Knepley 
210530ad9308SMatthew Knepley    Input Parameter:
210630ad9308SMatthew Knepley .  pc - the preconditioner context
210730ad9308SMatthew Knepley 
210830ad9308SMatthew Knepley    Output Parameters:
2109a04f6461SBarry Smith +  A00 - the (0,0) block
2110a04f6461SBarry Smith .  A01 - the (0,1) block
2111a04f6461SBarry Smith .  A10 - the (1,0) block
2112a04f6461SBarry Smith -  A11 - the (1,1) block
211330ad9308SMatthew Knepley 
211430ad9308SMatthew Knepley    Level: advanced
211530ad9308SMatthew Knepley 
211630ad9308SMatthew Knepley .seealso: PCFIELDSPLIT
211730ad9308SMatthew Knepley @*/
2118a04f6461SBarry Smith PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
211930ad9308SMatthew Knepley {
212030ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;
212130ad9308SMatthew Knepley 
212230ad9308SMatthew Knepley   PetscFunctionBegin;
21230700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2124ce94432eSBarry Smith   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
2125a04f6461SBarry Smith   if (A00) *A00 = jac->pmat[0];
2126a04f6461SBarry Smith   if (A01) *A01 = jac->B;
2127a04f6461SBarry Smith   if (A10) *A10 = jac->C;
2128a04f6461SBarry Smith   if (A11) *A11 = jac->pmat[1];
212930ad9308SMatthew Knepley   PetscFunctionReturn(0);
213030ad9308SMatthew Knepley }
213130ad9308SMatthew Knepley 
213279416396SBarry Smith #undef __FUNCT__
213379416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
21341e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
213579416396SBarry Smith {
213679416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
2137e69d4d44SBarry Smith   PetscErrorCode ierr;
213879416396SBarry Smith 
213979416396SBarry Smith   PetscFunctionBegin;
214079416396SBarry Smith   jac->type = type;
21413b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
21423b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
21433b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
21442fa5cd67SKarl Rupp 
2145bdf89e91SBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
214629f8a81cSJed Brown     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",PCFieldSplitSetSchurPre_FieldSplit);CHKERRQ(ierr);
214737a82bf0SJed Brown     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",PCFieldSplitGetSchurPre_FieldSplit);CHKERRQ(ierr);
2148bdf89e91SBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr);
2149e69d4d44SBarry Smith 
21503b224e63SBarry Smith   } else {
21513b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
21523b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
21532fa5cd67SKarl Rupp 
2154bdf89e91SBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
215529f8a81cSJed Brown     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurPre_C",0);CHKERRQ(ierr);
215637a82bf0SJed Brown     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSchurPre_C",0);CHKERRQ(ierr);
2157bdf89e91SBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",0);CHKERRQ(ierr);
21583b224e63SBarry Smith   }
215979416396SBarry Smith   PetscFunctionReturn(0);
216079416396SBarry Smith }
216179416396SBarry Smith 
216251f519a2SBarry Smith #undef __FUNCT__
216351f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
21641e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
216551f519a2SBarry Smith {
216651f519a2SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
216751f519a2SBarry Smith 
216851f519a2SBarry Smith   PetscFunctionBegin;
2169ce94432eSBarry Smith   if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
2170ce94432eSBarry 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);
217151f519a2SBarry Smith   jac->bs = bs;
217251f519a2SBarry Smith   PetscFunctionReturn(0);
217351f519a2SBarry Smith }
217451f519a2SBarry Smith 
217579416396SBarry Smith #undef __FUNCT__
217679416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
2177bc08b0f1SBarry Smith /*@
217879416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
217979416396SBarry Smith 
218079416396SBarry Smith    Collective on PC
218179416396SBarry Smith 
218279416396SBarry Smith    Input Parameter:
218379416396SBarry Smith .  pc - the preconditioner context
218481540f2fSBarry Smith .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
218579416396SBarry Smith 
218679416396SBarry Smith    Options Database Key:
2187a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
218879416396SBarry Smith 
2189b02e2d75SMatthew G Knepley    Level: Intermediate
219079416396SBarry Smith 
219179416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
219279416396SBarry Smith 
219379416396SBarry Smith .seealso: PCCompositeSetType()
219479416396SBarry Smith 
219579416396SBarry Smith @*/
21967087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
219779416396SBarry Smith {
21984ac538c5SBarry Smith   PetscErrorCode ierr;
219979416396SBarry Smith 
220079416396SBarry Smith   PetscFunctionBegin;
22010700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
22024ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr);
220379416396SBarry Smith   PetscFunctionReturn(0);
220479416396SBarry Smith }
220579416396SBarry Smith 
2206b02e2d75SMatthew G Knepley #undef __FUNCT__
2207b02e2d75SMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetType"
2208b02e2d75SMatthew G Knepley /*@
2209b02e2d75SMatthew G Knepley   PCFieldSplitGetType - Gets the type of fieldsplit preconditioner.
2210b02e2d75SMatthew G Knepley 
2211b02e2d75SMatthew G Knepley   Not collective
2212b02e2d75SMatthew G Knepley 
2213b02e2d75SMatthew G Knepley   Input Parameter:
2214b02e2d75SMatthew G Knepley . pc - the preconditioner context
2215b02e2d75SMatthew G Knepley 
2216b02e2d75SMatthew G Knepley   Output Parameter:
2217b02e2d75SMatthew G Knepley . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
2218b02e2d75SMatthew G Knepley 
2219b02e2d75SMatthew G Knepley   Level: Intermediate
2220b02e2d75SMatthew G Knepley 
2221b02e2d75SMatthew G Knepley .keywords: PC, set, type, composite preconditioner, additive, multiplicative
2222b02e2d75SMatthew G Knepley .seealso: PCCompositeSetType()
2223b02e2d75SMatthew G Knepley @*/
2224b02e2d75SMatthew G Knepley PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type)
2225b02e2d75SMatthew G Knepley {
2226b02e2d75SMatthew G Knepley   PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;
2227b02e2d75SMatthew G Knepley 
2228b02e2d75SMatthew G Knepley   PetscFunctionBegin;
2229b02e2d75SMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2230b02e2d75SMatthew G Knepley   PetscValidIntPointer(type,2);
2231b02e2d75SMatthew G Knepley   *type = jac->type;
2232b02e2d75SMatthew G Knepley   PetscFunctionReturn(0);
2233b02e2d75SMatthew G Knepley }
2234b02e2d75SMatthew G Knepley 
22354ab8060aSDmitry Karpeev #undef __FUNCT__
22364ab8060aSDmitry Karpeev #define __FUNCT__ "PCFieldSplitSetDMSplits"
22374ab8060aSDmitry Karpeev /*@
22384ab8060aSDmitry Karpeev    PCFieldSplitSetDMSplits - Flags whether DMCreateFieldDecomposition() should be used to define the splits, whenever possible.
22394ab8060aSDmitry Karpeev 
22404ab8060aSDmitry Karpeev    Logically Collective
22414ab8060aSDmitry Karpeev 
22424ab8060aSDmitry Karpeev    Input Parameters:
22434ab8060aSDmitry Karpeev +  pc   - the preconditioner context
22444ab8060aSDmitry Karpeev -  flg  - boolean indicating whether to use field splits defined by the DM
22454ab8060aSDmitry Karpeev 
22464ab8060aSDmitry Karpeev    Options Database Key:
22474ab8060aSDmitry Karpeev .  -pc_fieldsplit_dm_splits
22484ab8060aSDmitry Karpeev 
22494ab8060aSDmitry Karpeev    Level: Intermediate
22504ab8060aSDmitry Karpeev 
22514ab8060aSDmitry Karpeev .keywords: PC, DM, composite preconditioner, additive, multiplicative
22524ab8060aSDmitry Karpeev 
22534ab8060aSDmitry Karpeev .seealso: PCFieldSplitGetDMSplits()
22544ab8060aSDmitry Karpeev 
22554ab8060aSDmitry Karpeev @*/
22564ab8060aSDmitry Karpeev PetscErrorCode  PCFieldSplitSetDMSplits(PC pc,PetscBool flg)
22574ab8060aSDmitry Karpeev {
22584ab8060aSDmitry Karpeev   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
22594ab8060aSDmitry Karpeev   PetscBool      isfs;
22604ab8060aSDmitry Karpeev   PetscErrorCode ierr;
22614ab8060aSDmitry Karpeev 
22624ab8060aSDmitry Karpeev   PetscFunctionBegin;
22634ab8060aSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
22644ab8060aSDmitry Karpeev   PetscValidLogicalCollectiveBool(pc,flg,2);
22654ab8060aSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
22664ab8060aSDmitry Karpeev   if (isfs) {
22674ab8060aSDmitry Karpeev     jac->dm_splits = flg;
22684ab8060aSDmitry Karpeev   }
22694ab8060aSDmitry Karpeev   PetscFunctionReturn(0);
22704ab8060aSDmitry Karpeev }
22714ab8060aSDmitry Karpeev 
22724ab8060aSDmitry Karpeev 
22734ab8060aSDmitry Karpeev #undef __FUNCT__
22744ab8060aSDmitry Karpeev #define __FUNCT__ "PCFieldSplitGetDMSplits"
22754ab8060aSDmitry Karpeev /*@
22764ab8060aSDmitry Karpeev    PCFieldSplitGetDMSplits - Returns flag indicating whether DMCreateFieldDecomposition() should be used to define the splits, whenever possible.
22774ab8060aSDmitry Karpeev 
22784ab8060aSDmitry Karpeev    Logically Collective
22794ab8060aSDmitry Karpeev 
22804ab8060aSDmitry Karpeev    Input Parameter:
22814ab8060aSDmitry Karpeev .  pc   - the preconditioner context
22824ab8060aSDmitry Karpeev 
22834ab8060aSDmitry Karpeev    Output Parameter:
22844ab8060aSDmitry Karpeev .  flg  - boolean indicating whether to use field splits defined by the DM
22854ab8060aSDmitry Karpeev 
22864ab8060aSDmitry Karpeev    Level: Intermediate
22874ab8060aSDmitry Karpeev 
22884ab8060aSDmitry Karpeev .keywords: PC, DM, composite preconditioner, additive, multiplicative
22894ab8060aSDmitry Karpeev 
22904ab8060aSDmitry Karpeev .seealso: PCFieldSplitSetDMSplits()
22914ab8060aSDmitry Karpeev 
22924ab8060aSDmitry Karpeev @*/
22934ab8060aSDmitry Karpeev PetscErrorCode  PCFieldSplitGetDMSplits(PC pc,PetscBool* flg)
22944ab8060aSDmitry Karpeev {
22954ab8060aSDmitry Karpeev   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
22964ab8060aSDmitry Karpeev   PetscBool      isfs;
22974ab8060aSDmitry Karpeev   PetscErrorCode ierr;
22984ab8060aSDmitry Karpeev 
22994ab8060aSDmitry Karpeev   PetscFunctionBegin;
23004ab8060aSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
23014ab8060aSDmitry Karpeev   PetscValidPointer(flg,2);
23024ab8060aSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
23034ab8060aSDmitry Karpeev   if (isfs) {
23044ab8060aSDmitry Karpeev     if(flg) *flg = jac->dm_splits;
23054ab8060aSDmitry Karpeev   }
23064ab8060aSDmitry Karpeev   PetscFunctionReturn(0);
23074ab8060aSDmitry Karpeev }
23084ab8060aSDmitry Karpeev 
23090971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
23100971522cSBarry Smith /*MC
2311a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
2312a04f6461SBarry Smith                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.
23130971522cSBarry Smith 
2314edf189efSBarry Smith      To set options on the solvers for each block append -fieldsplit_ to all the PC
2315edf189efSBarry Smith         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
23160971522cSBarry Smith 
2317a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
231869a612a9SBarry Smith          and set the options directly on the resulting KSP object
23190971522cSBarry Smith 
23200971522cSBarry Smith    Level: intermediate
23210971522cSBarry Smith 
232279416396SBarry Smith    Options Database Keys:
232381540f2fSBarry Smith +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
232481540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
232581540f2fSBarry Smith                               been supplied explicitly by -pc_fieldsplit_%d_fields
232681540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
23270f188ba9SJed Brown .   -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur> - type of relaxation or factorization splitting
2328a6a584a2SBarry Smith .   -pc_fieldsplit_schur_precondition <self,selfp,user,a11,full> - default is a11; see PCFieldSplitSetSchurPre()
2329435f959eSBarry Smith .   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver
233079416396SBarry Smith 
23315d4c12cdSJungho Lee -    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
23325d4c12cdSJungho Lee      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
23335d4c12cdSJungho Lee 
2334c8a0d604SMatthew G Knepley    Notes:
2335c8a0d604SMatthew G Knepley     Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
2336d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
2337d32f9abdSBarry Smith 
2338d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
2339d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
2340d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
2341d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
2342d32f9abdSBarry Smith 
2343c8a0d604SMatthew G Knepley $     For the Schur complement preconditioner if J = ( A00 A01 )
2344c8a0d604SMatthew G Knepley $                                                    ( A10 A11 )
2345c8a0d604SMatthew G Knepley $     the preconditioner using full factorization is
234613b05affSJed Brown $              ( I   -ksp(A00) A01 ) ( inv(A00)     0  ) (     I          0  )
2347c8a0d604SMatthew G Knepley $              ( 0         I       ) (   0      ksp(S) ) ( -A10 ksp(A00)  I  )
234813b05affSJed Brown      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_.  S is the Schur complement
234913b05affSJed Brown $              S = A11 - A10 ksp(A00) A01
235013b05affSJed 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
2351c8a0d604SMatthew G Knepley      in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0,
2352c8a0d604SMatthew G Knepley      it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
2353a6a584a2SBarry Smith      A11 is used to construct a preconditioner for S, use PCFieldSplitSetSchurPre() for all the possible ways to construct the preconditioner for S.
23541d71051fSDmitry Karpeev 
2355a7476a74SDmitry Karpeev      The factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above,
23565668aaf4SBarry Smith      diag gives
2357c8a0d604SMatthew G Knepley $              ( inv(A00)     0   )
2358c8a0d604SMatthew G Knepley $              (   0      -ksp(S) )
23595668aaf4SBarry 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
2360c8a0d604SMatthew G Knepley $              (  A00   0 )
2361c8a0d604SMatthew G Knepley $              (  A10   S )
2362c8a0d604SMatthew G Knepley      where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of
2363c8a0d604SMatthew G Knepley $              ( A00 A01 )
2364c8a0d604SMatthew G Knepley $              (  0   S  )
2365c8a0d604SMatthew G Knepley      where again the inverses of A00 and S are applied using KSPs.
2366e69d4d44SBarry Smith 
2367edf189efSBarry Smith      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
2368edf189efSBarry Smith      is used automatically for a second block.
2369edf189efSBarry Smith 
2370ff218e97SBarry Smith      The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1.
2371ff218e97SBarry Smith      Generally it should be used with the AIJ format.
2372ff218e97SBarry Smith 
2373ff218e97SBarry Smith      The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see,
2374ff218e97SBarry Smith      for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT
2375ff218e97SBarry Smith      inside a smoother resulting in "Distributive Smoothers".
23760716a85fSBarry Smith 
2377a541d17aSBarry Smith    Concepts: physics based preconditioners, block preconditioners
23780971522cSBarry Smith 
23793bc631e0SBarry Smith    There is a nice discussion of block preconditioners in
23803bc631e0SBarry Smith 
23813bc631e0SBarry Smith [El08] A taxonomy and comparison of parallel block multi-level preconditioners for the incompressible Navier-Stokes equations
23823bc631e0SBarry Smith        Howard Elman, V.E. Howle, John Shadid, Robert Shuttleworth, Ray Tuminaro, Journal of Computational Physics 227 (2008) 1790--1808
23833bc631e0SBarry Smith        http://chess.cs.umd.edu/~elman/papers/tax.pdf
23843bc631e0SBarry Smith 
2385a6a584a2SBarry 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
2386a6a584a2SBarry 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
2387a6a584a2SBarry Smith    of backsolving for the saturations in the last step it solves a full coupled (ILU) system for updates to all the variables.
2388a6a584a2SBarry Smith 
23897e8cb189SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
23901d71051fSDmitry Karpeev            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSetSchurPre(),
23911d71051fSDmitry Karpeev 	   MatSchurComplementSetAinvType()
23920971522cSBarry Smith M*/
23930971522cSBarry Smith 
23940971522cSBarry Smith #undef __FUNCT__
23950971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
23968cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_FieldSplit(PC pc)
23970971522cSBarry Smith {
23980971522cSBarry Smith   PetscErrorCode ierr;
23990971522cSBarry Smith   PC_FieldSplit  *jac;
24000971522cSBarry Smith 
24010971522cSBarry Smith   PetscFunctionBegin;
2402b00a9115SJed Brown   ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr);
24032fa5cd67SKarl Rupp 
24040971522cSBarry Smith   jac->bs                 = -1;
24050971522cSBarry Smith   jac->nsplits            = 0;
24063e197d65SBarry Smith   jac->type               = PC_COMPOSITE_MULTIPLICATIVE;
2407e6cab6aaSJed Brown   jac->schurpre           = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
2408c9c6ffaaSJed Brown   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL;
2409fbe7908bSJed Brown   jac->dm_splits          = PETSC_TRUE;
241051f519a2SBarry Smith 
24110971522cSBarry Smith   pc->data = (void*)jac;
24120971522cSBarry Smith 
24130971522cSBarry Smith   pc->ops->apply           = PCApply_FieldSplit;
2414421e10b8SBarry Smith   pc->ops->applytranspose  = PCApplyTranspose_FieldSplit;
24150971522cSBarry Smith   pc->ops->setup           = PCSetUp_FieldSplit;
2416574deadeSBarry Smith   pc->ops->reset           = PCReset_FieldSplit;
24170971522cSBarry Smith   pc->ops->destroy         = PCDestroy_FieldSplit;
24180971522cSBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_FieldSplit;
24190971522cSBarry Smith   pc->ops->view            = PCView_FieldSplit;
24200971522cSBarry Smith   pc->ops->applyrichardson = 0;
24210971522cSBarry Smith 
2422bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
2423bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
2424bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
2425bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
2426bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
24276dbb499eSCian Wilson   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitRestrictIS_C",PCFieldSplitRestrictIS_FieldSplit);CHKERRQ(ierr);
24280971522cSBarry Smith   PetscFunctionReturn(0);
24290971522cSBarry Smith }
24300971522cSBarry Smith 
24310971522cSBarry Smith 
2432a541d17aSBarry Smith 
2433