xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 2e71c61d18e71e694bfbe2b314256bf3aa85754d)
1dba47a55SKris Buschelman 
2b45d2f2cSJed Brown #include <petsc-private/pcimpl.h>     /*I "petscpc.h" I*/
31e25c274SJed Brown #include <petscdm.h>
40971522cSBarry Smith 
5443836d0SMatthew G Knepley /*
6443836d0SMatthew G Knepley   There is a nice discussion of block preconditioners in
7443836d0SMatthew G Knepley 
81997fe2eSSatish Balay [El08] A taxonomy and comparison of parallel block multi-level preconditioners for the incompressible Navier-Stokes equations
9443836d0SMatthew G Knepley        Howard Elman, V.E. Howle, John Shadid, Robert Shuttleworth, Ray Tuminaro, Journal of Computational Physics 227 (2008) 1790--1808
101b8e4c5fSJed Brown        http://chess.cs.umd.edu/~elman/papers/tax.pdf
11443836d0SMatthew G Knepley */
12443836d0SMatthew G Knepley 
13*2e71c61dSMatthew G. Knepley const char *const PCFieldSplitSchurPreTypes[] = {"SELF","SELFP","A11","USER","FULL","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0};
14c9c6ffaaSJed Brown const char *const PCFieldSplitSchurFactTypes[] = {"DIAG","LOWER","UPPER","FULL","PCFieldSplitSchurFactType","PC_FIELDSPLIT_SCHUR_FACT_",0};
15c5d2311dSJed Brown 
160971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
170971522cSBarry Smith struct _PC_FieldSplitLink {
1869a612a9SBarry Smith   KSP               ksp;
19443836d0SMatthew G Knepley   Vec               x,y,z;
20db4c96c1SJed Brown   char              *splitname;
210971522cSBarry Smith   PetscInt          nfields;
225d4c12cdSJungho Lee   PetscInt          *fields,*fields_col;
231b9fc7fcSBarry Smith   VecScatter        sctx;
245d4c12cdSJungho Lee   IS                is,is_col;
2551f519a2SBarry Smith   PC_FieldSplitLink next,previous;
260971522cSBarry Smith };
270971522cSBarry Smith 
280971522cSBarry Smith typedef struct {
2968dd23aaSBarry Smith   PCCompositeType type;
30ace3abfcSBarry Smith   PetscBool       defaultsplit;                    /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */
31ace3abfcSBarry Smith   PetscBool       splitdefined;                    /* Flag is set after the splits have been defined, to prevent more splits from being added */
3230ad9308SMatthew Knepley   PetscInt        bs;                              /* Block size for IS and Mat structures */
3330ad9308SMatthew Knepley   PetscInt        nsplits;                         /* Number of field divisions defined */
3479416396SBarry Smith   Vec             *x,*y,w1,w2;
35519d70e2SJed Brown   Mat             *mat;                            /* The diagonal block for each split */
36519d70e2SJed Brown   Mat             *pmat;                           /* The preconditioning diagonal block for each split */
3730ad9308SMatthew Knepley   Mat             *Afield;                         /* The rows of the matrix associated with each split */
38ace3abfcSBarry Smith   PetscBool       issetup;
392fa5cd67SKarl Rupp 
4030ad9308SMatthew Knepley   /* Only used when Schur complement preconditioning is used */
4130ad9308SMatthew Knepley   Mat                       B;                     /* The (0,1) block */
4230ad9308SMatthew Knepley   Mat                       C;                     /* The (1,0) block */
43443836d0SMatthew G Knepley   Mat                       schur;                 /* The Schur complement S = A11 - A10 A00^{-1} A01, the KSP here, kspinner, is H_1 in [El08] */
44a7476a74SDmitry Karpeev   Mat                       schurp;                /* Assembled approximation to S built by MatSchurComplement to be used as a preconditioning matrix when solving with S */
45084e4875SJed Brown   Mat                       schur_user;            /* User-provided preconditioning matrix for the Schur complement */
46084e4875SJed Brown   PCFieldSplitSchurPreType  schurpre;              /* Determines which preconditioning matrix is used for the Schur complement */
47c9c6ffaaSJed Brown   PCFieldSplitSchurFactType schurfactorization;
4830ad9308SMatthew Knepley   KSP                       kspschur;              /* The solver for S */
49443836d0SMatthew G Knepley   KSP                       kspupper;              /* The solver for A in the upper diagonal part of the factorization (H_2 in [El08]) */
5097bbdb24SBarry Smith   PC_FieldSplitLink         head;
5163ec74ffSBarry Smith   PetscBool                 reset;                  /* indicates PCReset() has been last called on this object, hack */
52c1570756SJed Brown   PetscBool                 suboptionsset;          /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */
534ab8060aSDmitry Karpeev   PetscBool                 dm_splits;              /* Whether to use DMCreateFieldDecomposition() whenever possible */
540971522cSBarry Smith } PC_FieldSplit;
550971522cSBarry Smith 
5616913363SBarry Smith /*
5716913363SBarry Smith     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
5816913363SBarry Smith    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
5916913363SBarry Smith    PC you could change this.
6016913363SBarry Smith */
61084e4875SJed Brown 
62e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it.  This way the
63084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */
64084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac)
65084e4875SJed Brown {
66084e4875SJed Brown   switch (jac->schurpre) {
67084e4875SJed Brown   case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur;
68a7476a74SDmitry Karpeev   case PC_FIELDSPLIT_SCHUR_PRE_SELFP: return jac->schurp;
69e87fab1bSBarry Smith   case PC_FIELDSPLIT_SCHUR_PRE_A11: return jac->pmat[1];
70e74569cdSMatthew G. Knepley   case PC_FIELDSPLIT_SCHUR_PRE_FULL: /* We calculate this and store it in schur_user */
71084e4875SJed Brown   case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */
72084e4875SJed Brown   default:
73084e4875SJed Brown     return jac->schur_user ? jac->schur_user : jac->pmat[1];
74084e4875SJed Brown   }
75084e4875SJed Brown }
76084e4875SJed Brown 
77084e4875SJed Brown 
789804daf3SBarry Smith #include <petscdraw.h>
790971522cSBarry Smith #undef __FUNCT__
800971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
810971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
820971522cSBarry Smith {
830971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
840971522cSBarry Smith   PetscErrorCode    ierr;
85d9884438SBarry Smith   PetscBool         iascii,isdraw;
860971522cSBarry Smith   PetscInt          i,j;
875a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
880971522cSBarry Smith 
890971522cSBarry Smith   PetscFunctionBegin;
90251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
91d9884438SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
920971522cSBarry Smith   if (iascii) {
932c9966d7SBarry Smith     if (jac->bs > 0) {
9451f519a2SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
952c9966d7SBarry Smith     } else {
962c9966d7SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D\n",PCCompositeTypes[jac->type],jac->nsplits);CHKERRQ(ierr);
972c9966d7SBarry Smith     }
98f5236f50SJed Brown     if (pc->useAmat) {
99f5236f50SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for blocks\n");CHKERRQ(ierr);
100a3df900dSMatthew G Knepley     }
10169a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
1020971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1030971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
1041ab39975SBarry Smith       if (ilink->fields) {
1050971522cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
10679416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
1075a9f2f41SSatish Balay         for (j=0; j<ilink->nfields; j++) {
10879416396SBarry Smith           if (j > 0) {
10979416396SBarry Smith             ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
11079416396SBarry Smith           }
1115a9f2f41SSatish Balay           ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
1120971522cSBarry Smith         }
1130971522cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
11479416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1151ab39975SBarry Smith       } else {
1161ab39975SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1171ab39975SBarry Smith       }
1185a9f2f41SSatish Balay       ierr  = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
1195a9f2f41SSatish Balay       ilink = ilink->next;
1200971522cSBarry Smith     }
1210971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1222fa5cd67SKarl Rupp   }
1232fa5cd67SKarl Rupp 
1242fa5cd67SKarl Rupp  if (isdraw) {
125d9884438SBarry Smith     PetscDraw draw;
126d9884438SBarry Smith     PetscReal x,y,w,wd;
127d9884438SBarry Smith 
128d9884438SBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
129d9884438SBarry Smith     ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
130d9884438SBarry Smith     w    = 2*PetscMin(1.0 - x,x);
131d9884438SBarry Smith     wd   = w/(jac->nsplits + 1);
132d9884438SBarry Smith     x    = x - wd*(jac->nsplits-1)/2.0;
133d9884438SBarry Smith     for (i=0; i<jac->nsplits; i++) {
134d9884438SBarry Smith       ierr  = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
135d9884438SBarry Smith       ierr  = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
136d9884438SBarry Smith       ierr  = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
137d9884438SBarry Smith       x    += wd;
138d9884438SBarry Smith       ilink = ilink->next;
139d9884438SBarry Smith     }
1400971522cSBarry Smith   }
1410971522cSBarry Smith   PetscFunctionReturn(0);
1420971522cSBarry Smith }
1430971522cSBarry Smith 
1440971522cSBarry Smith #undef __FUNCT__
1453b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur"
1463b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
1473b224e63SBarry Smith {
1483b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1493b224e63SBarry Smith   PetscErrorCode    ierr;
1504996c5bdSBarry Smith   PetscBool         iascii,isdraw;
1513b224e63SBarry Smith   PetscInt          i,j;
1523b224e63SBarry Smith   PC_FieldSplitLink ilink = jac->head;
1533b224e63SBarry Smith 
1543b224e63SBarry Smith   PetscFunctionBegin;
155251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1564996c5bdSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
1573b224e63SBarry Smith   if (iascii) {
1582c9966d7SBarry Smith     if (jac->bs > 0) {
159c9c6ffaaSJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
1602c9966d7SBarry Smith     } else {
1612c9966d7SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, factorization %s\n",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
1622c9966d7SBarry Smith     }
163f5236f50SJed Brown     if (pc->useAmat) {
164f5236f50SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  using Amat (not Pmat) as operator for blocks\n");CHKERRQ(ierr);
165a3df900dSMatthew G Knepley     }
1663e8b8b31SMatthew G Knepley     switch (jac->schurpre) {
1673e8b8b31SMatthew G Knepley     case PC_FIELDSPLIT_SCHUR_PRE_SELF:
1683e8b8b31SMatthew G Knepley       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from S itself\n");CHKERRQ(ierr);break;
169a7476a74SDmitry Karpeev     case PC_FIELDSPLIT_SCHUR_PRE_SELFP:
170a7476a74SDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from Sp, an assembled approximation to S, which uses A00's diagonal's inverse\n");CHKERRQ(ierr);break;
171e87fab1bSBarry Smith     case PC_FIELDSPLIT_SCHUR_PRE_A11:
172e87fab1bSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from A11\n");CHKERRQ(ierr);break;
173e74569cdSMatthew G. Knepley     case PC_FIELDSPLIT_SCHUR_PRE_FULL:
174e74569cdSMatthew G. Knepley       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from the exact Schur complement\n");CHKERRQ(ierr);break;
1753e8b8b31SMatthew G Knepley     case PC_FIELDSPLIT_SCHUR_PRE_USER:
1763e8b8b31SMatthew G Knepley       if (jac->schur_user) {
1773e8b8b31SMatthew G Knepley         ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from user provided matrix\n");CHKERRQ(ierr);
1783e8b8b31SMatthew G Knepley       } else {
179e87fab1bSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from A11\n");CHKERRQ(ierr);
1803e8b8b31SMatthew G Knepley       }
1813e8b8b31SMatthew G Knepley       break;
1823e8b8b31SMatthew G Knepley     default:
18382f516ccSBarry Smith       SETERRQ1(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Invalid Schur preconditioning type: %d", jac->schurpre);
1843e8b8b31SMatthew G Knepley     }
1853b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Split info:\n");CHKERRQ(ierr);
1863b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1873b224e63SBarry Smith     for (i=0; i<jac->nsplits; i++) {
1883b224e63SBarry Smith       if (ilink->fields) {
1893b224e63SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
1903b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
1913b224e63SBarry Smith         for (j=0; j<ilink->nfields; j++) {
1923b224e63SBarry Smith           if (j > 0) {
1933b224e63SBarry Smith             ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
1943b224e63SBarry Smith           }
1953b224e63SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
1963b224e63SBarry Smith         }
1973b224e63SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1983b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1993b224e63SBarry Smith       } else {
2003b224e63SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
2013b224e63SBarry Smith       }
2023b224e63SBarry Smith       ilink = ilink->next;
2033b224e63SBarry Smith     }
204435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block\n");CHKERRQ(ierr);
2053b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
20606de4afeSJed Brown     if (jac->head) {
207443836d0SMatthew G Knepley       ierr = KSPView(jac->head->ksp,viewer);CHKERRQ(ierr);
20806de4afeSJed Brown     } else  {ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);}
2093b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
21006de4afeSJed Brown     if (jac->head && jac->kspupper != jac->head->ksp) {
211443836d0SMatthew G Knepley       ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for upper A00 in upper triangular factor \n");CHKERRQ(ierr);
212443836d0SMatthew G Knepley       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
213b2750c55SJed Brown       if (jac->kspupper) {ierr = KSPView(jac->kspupper,viewer);CHKERRQ(ierr);}
214b2750c55SJed Brown       else {ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);}
215443836d0SMatthew G Knepley       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
216443836d0SMatthew G Knepley     }
217435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr);
2183b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
21912cae6f2SJed Brown     if (jac->kspschur) {
2203b224e63SBarry Smith       ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
22112cae6f2SJed Brown     } else {
22212cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
22312cae6f2SJed Brown     }
2243b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
2253b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
22606de4afeSJed Brown   } else if (isdraw && jac->head) {
2274996c5bdSBarry Smith     PetscDraw draw;
2284996c5bdSBarry Smith     PetscReal x,y,w,wd,h;
2294996c5bdSBarry Smith     PetscInt  cnt = 2;
2304996c5bdSBarry Smith     char      str[32];
2314996c5bdSBarry Smith 
2324996c5bdSBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
2334996c5bdSBarry Smith     ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
234c74581afSBarry Smith     if (jac->kspupper != jac->head->ksp) cnt++;
235c74581afSBarry Smith     w  = 2*PetscMin(1.0 - x,x);
236c74581afSBarry Smith     wd = w/(cnt + 1);
237c74581afSBarry Smith 
2384996c5bdSBarry Smith     ierr = PetscSNPrintf(str,32,"Schur fact. %s",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
2390298fd71SBarry Smith     ierr = PetscDrawBoxedString(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr);
2404996c5bdSBarry Smith     y   -= h;
2414996c5bdSBarry Smith     if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_USER &&  !jac->schur_user) {
242e87fab1bSBarry Smith       ierr = PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[PC_FIELDSPLIT_SCHUR_PRE_A11]);CHKERRQ(ierr);
2433b224e63SBarry Smith     } else {
2444996c5bdSBarry Smith       ierr = PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[jac->schurpre]);CHKERRQ(ierr);
2454996c5bdSBarry Smith     }
2460298fd71SBarry Smith     ierr = PetscDrawBoxedString(draw,x+wd*(cnt-1)/2.0,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr);
2474996c5bdSBarry Smith     y   -= h;
2484996c5bdSBarry Smith     x    = x - wd*(cnt-1)/2.0;
2494996c5bdSBarry Smith 
2504996c5bdSBarry Smith     ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
2514996c5bdSBarry Smith     ierr = KSPView(jac->head->ksp,viewer);CHKERRQ(ierr);
2524996c5bdSBarry Smith     ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
2534996c5bdSBarry Smith     if (jac->kspupper != jac->head->ksp) {
2544996c5bdSBarry Smith       x   += wd;
2554996c5bdSBarry Smith       ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
2564996c5bdSBarry Smith       ierr = KSPView(jac->kspupper,viewer);CHKERRQ(ierr);
2574996c5bdSBarry Smith       ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
2584996c5bdSBarry Smith     }
2594996c5bdSBarry Smith     x   += wd;
2604996c5bdSBarry Smith     ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
2614996c5bdSBarry Smith     ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
2624996c5bdSBarry Smith     ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
2633b224e63SBarry Smith   }
2643b224e63SBarry Smith   PetscFunctionReturn(0);
2653b224e63SBarry Smith }
2663b224e63SBarry Smith 
2673b224e63SBarry Smith #undef __FUNCT__
2686c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private"
2696c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */
2706c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc)
2716c924f48SJed Brown {
2726c924f48SJed Brown   PetscErrorCode ierr;
2736c924f48SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
2745d4c12cdSJungho Lee   PetscInt       i,nfields,*ifields,nfields_col,*ifields_col;
2755d4c12cdSJungho Lee   PetscBool      flg,flg_col;
2765d4c12cdSJungho Lee   char           optionname[128],splitname[8],optionname_col[128];
2776c924f48SJed Brown 
2786c924f48SJed Brown   PetscFunctionBegin;
279785e854fSJed Brown   ierr = PetscMalloc1(jac->bs,&ifields);CHKERRQ(ierr);
280785e854fSJed Brown   ierr = PetscMalloc1(jac->bs,&ifields_col);CHKERRQ(ierr);
2816c924f48SJed Brown   for (i=0,flg=PETSC_TRUE;; i++) {
2828caf3d72SBarry Smith     ierr        = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr);
2838caf3d72SBarry Smith     ierr        = PetscSNPrintf(optionname,sizeof(optionname),"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr);
2848caf3d72SBarry Smith     ierr        = PetscSNPrintf(optionname_col,sizeof(optionname_col),"-pc_fieldsplit_%D_fields_col",i);CHKERRQ(ierr);
2856c924f48SJed Brown     nfields     = jac->bs;
28629499fbbSJungho Lee     nfields_col = jac->bs;
2876c924f48SJed Brown     ierr        = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr);
2885d4c12cdSJungho Lee     ierr        = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);CHKERRQ(ierr);
2896c924f48SJed Brown     if (!flg) break;
2905d4c12cdSJungho Lee     else if (flg && !flg_col) {
2915d4c12cdSJungho Lee       if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
2925d4c12cdSJungho Lee       ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);CHKERRQ(ierr);
2932fa5cd67SKarl Rupp     } else {
2945d4c12cdSJungho Lee       if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
2955d4c12cdSJungho Lee       if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match");
2965d4c12cdSJungho Lee       ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);CHKERRQ(ierr);
2975d4c12cdSJungho Lee     }
2986c924f48SJed Brown   }
2996c924f48SJed Brown   if (i > 0) {
3006c924f48SJed Brown     /* Makes command-line setting of splits take precedence over setting them in code.
3016c924f48SJed Brown        Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would
3026c924f48SJed Brown        create new splits, which would probably not be what the user wanted. */
3036c924f48SJed Brown     jac->splitdefined = PETSC_TRUE;
3046c924f48SJed Brown   }
3056c924f48SJed Brown   ierr = PetscFree(ifields);CHKERRQ(ierr);
3065d4c12cdSJungho Lee   ierr = PetscFree(ifields_col);CHKERRQ(ierr);
3076c924f48SJed Brown   PetscFunctionReturn(0);
3086c924f48SJed Brown }
3096c924f48SJed Brown 
3106c924f48SJed Brown #undef __FUNCT__
31169a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
31269a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
3130971522cSBarry Smith {
3140971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
3150971522cSBarry Smith   PetscErrorCode    ierr;
3165a9f2f41SSatish Balay   PC_FieldSplitLink ilink              = jac->head;
3177287d2a3SDmitry Karpeev   PetscBool         fieldsplit_default = PETSC_FALSE,stokes = PETSC_FALSE;
3186c924f48SJed Brown   PetscInt          i;
3190971522cSBarry Smith 
3200971522cSBarry Smith   PetscFunctionBegin;
3217287d2a3SDmitry Karpeev   /*
3227287d2a3SDmitry Karpeev    Kinda messy, but at least this now uses DMCreateFieldDecomposition() even with jac->reset.
3237287d2a3SDmitry Karpeev    Should probably be rewritten.
3247287d2a3SDmitry Karpeev    */
3257287d2a3SDmitry Karpeev   if (!ilink || jac->reset) {
3260298fd71SBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);CHKERRQ(ierr);
3274ab8060aSDmitry Karpeev     if (pc->dm && jac->dm_splits && !stokes) {
328bafc1b83SMatthew G Knepley       PetscInt  numFields, f, i, j;
3290784a22cSJed Brown       char      **fieldNames;
3307b62db95SJungho Lee       IS        *fields;
331e7c4fc90SDmitry Karpeev       DM        *dms;
332bafc1b83SMatthew G Knepley       DM        subdm[128];
333bafc1b83SMatthew G Knepley       PetscBool flg;
334bafc1b83SMatthew G Knepley 
33516621825SDmitry Karpeev       ierr = DMCreateFieldDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms);CHKERRQ(ierr);
336bafc1b83SMatthew G Knepley       /* Allow the user to prescribe the splits */
337bafc1b83SMatthew G Knepley       for (i = 0, flg = PETSC_TRUE;; i++) {
338bafc1b83SMatthew G Knepley         PetscInt ifields[128];
339bafc1b83SMatthew G Knepley         IS       compField;
340bafc1b83SMatthew G Knepley         char     optionname[128], splitname[8];
341bafc1b83SMatthew G Knepley         PetscInt nfields = numFields;
342bafc1b83SMatthew G Knepley 
3438caf3d72SBarry Smith         ierr = PetscSNPrintf(optionname, sizeof(optionname), "-pc_fieldsplit_%D_fields", i);CHKERRQ(ierr);
344bafc1b83SMatthew G Knepley         ierr = PetscOptionsGetIntArray(((PetscObject) pc)->prefix, optionname, ifields, &nfields, &flg);CHKERRQ(ierr);
345bafc1b83SMatthew G Knepley         if (!flg) break;
34682f516ccSBarry Smith         if (numFields > 128) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot currently support %d > 128 fields", numFields);
347bafc1b83SMatthew G Knepley         ierr = DMCreateSubDM(pc->dm, nfields, ifields, &compField, &subdm[i]);CHKERRQ(ierr);
348bafc1b83SMatthew G Knepley         if (nfields == 1) {
349bafc1b83SMatthew G Knepley           ierr = PCFieldSplitSetIS(pc, fieldNames[ifields[0]], compField);CHKERRQ(ierr);
35082f516ccSBarry Smith           /* ierr = PetscPrintf(PetscObjectComm((PetscObject)pc), "%s Field Indices:", fieldNames[ifields[0]]);CHKERRQ(ierr);
3510298fd71SBarry Smith              ierr = ISView(compField, NULL);CHKERRQ(ierr); */
352bafc1b83SMatthew G Knepley         } else {
3538caf3d72SBarry Smith           ierr = PetscSNPrintf(splitname, sizeof(splitname), "%D", i);CHKERRQ(ierr);
354bafc1b83SMatthew G Knepley           ierr = PCFieldSplitSetIS(pc, splitname, compField);CHKERRQ(ierr);
35582f516ccSBarry Smith           /* ierr = PetscPrintf(PetscObjectComm((PetscObject)pc), "%s Field Indices:", splitname);CHKERRQ(ierr);
3560298fd71SBarry Smith              ierr = ISView(compField, NULL);CHKERRQ(ierr); */
3577287d2a3SDmitry Karpeev         }
358bafc1b83SMatthew G Knepley         ierr = ISDestroy(&compField);CHKERRQ(ierr);
359bafc1b83SMatthew G Knepley         for (j = 0; j < nfields; ++j) {
360bafc1b83SMatthew G Knepley           f    = ifields[j];
3617b62db95SJungho Lee           ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr);
3627b62db95SJungho Lee           ierr = ISDestroy(&fields[f]);CHKERRQ(ierr);
3637b62db95SJungho Lee         }
364bafc1b83SMatthew G Knepley       }
365bafc1b83SMatthew G Knepley       if (i == 0) {
366bafc1b83SMatthew G Knepley         for (f = 0; f < numFields; ++f) {
367bafc1b83SMatthew G Knepley           ierr = PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);CHKERRQ(ierr);
368bafc1b83SMatthew G Knepley           ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr);
369bafc1b83SMatthew G Knepley           ierr = ISDestroy(&fields[f]);CHKERRQ(ierr);
370bafc1b83SMatthew G Knepley         }
371bafc1b83SMatthew G Knepley       } else {
372d724dfffSBarry Smith         for (j=0; j<numFields; j++) {
373d724dfffSBarry Smith           ierr = DMDestroy(dms+j);CHKERRQ(ierr);
374d724dfffSBarry Smith         }
375d724dfffSBarry Smith         ierr = PetscFree(dms);CHKERRQ(ierr);
376785e854fSJed Brown         ierr = PetscMalloc1(i, &dms);CHKERRQ(ierr);
3772fa5cd67SKarl Rupp         for (j = 0; j < i; ++j) dms[j] = subdm[j];
378bafc1b83SMatthew G Knepley       }
3797b62db95SJungho Lee       ierr = PetscFree(fieldNames);CHKERRQ(ierr);
3807b62db95SJungho Lee       ierr = PetscFree(fields);CHKERRQ(ierr);
381e7c4fc90SDmitry Karpeev       if (dms) {
3828b8307b2SJed Brown         ierr = PetscInfo(pc, "Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
383bafc1b83SMatthew G Knepley         for (ilink = jac->head, i = 0; ilink; ilink = ilink->next, ++i) {
3847287d2a3SDmitry Karpeev           const char *prefix;
3857287d2a3SDmitry Karpeev           ierr = PetscObjectGetOptionsPrefix((PetscObject)(ilink->ksp),&prefix);CHKERRQ(ierr);
3867287d2a3SDmitry Karpeev           ierr = PetscObjectSetOptionsPrefix((PetscObject)(dms[i]), prefix);CHKERRQ(ierr);
3877b62db95SJungho Lee           ierr = KSPSetDM(ilink->ksp, dms[i]);CHKERRQ(ierr);
3887b62db95SJungho Lee           ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE);CHKERRQ(ierr);
3897287d2a3SDmitry Karpeev           ierr = PetscObjectIncrementTabLevel((PetscObject)dms[i],(PetscObject)ilink->ksp,0);CHKERRQ(ierr);
390e7c4fc90SDmitry Karpeev           ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);
3912fa5ba8aSJed Brown         }
3927b62db95SJungho Lee         ierr = PetscFree(dms);CHKERRQ(ierr);
3938b8307b2SJed Brown       }
39466ffff09SJed Brown     } else {
395521d7252SBarry Smith       if (jac->bs <= 0) {
396704ba839SBarry Smith         if (pc->pmat) {
397521d7252SBarry Smith           ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
3982fa5cd67SKarl Rupp         } else jac->bs = 1;
399521d7252SBarry Smith       }
400d32f9abdSBarry Smith 
4016ce1633cSBarry Smith       if (stokes) {
4026ce1633cSBarry Smith         IS       zerodiags,rest;
4036ce1633cSBarry Smith         PetscInt nmin,nmax;
4046ce1633cSBarry Smith 
4056ce1633cSBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
4066ce1633cSBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
4076ce1633cSBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
4087287d2a3SDmitry Karpeev         if (jac->reset) {
4097287d2a3SDmitry Karpeev           jac->head->is       = rest;
4107287d2a3SDmitry Karpeev           jac->head->next->is = zerodiags;
4112fa5cd67SKarl Rupp         } else {
4126ce1633cSBarry Smith           ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
4136ce1633cSBarry Smith           ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr);
4147287d2a3SDmitry Karpeev         }
415fcfd50ebSBarry Smith         ierr = ISDestroy(&zerodiags);CHKERRQ(ierr);
416fcfd50ebSBarry Smith         ierr = ISDestroy(&rest);CHKERRQ(ierr);
4176ce1633cSBarry Smith       } else {
418ce94432eSBarry Smith         if (jac->reset) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used");
4199eeaaa73SJed Brown         ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&fieldsplit_default,NULL);CHKERRQ(ierr);
4207287d2a3SDmitry Karpeev         if (!fieldsplit_default) {
421d32f9abdSBarry Smith           /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
422d32f9abdSBarry Smith            then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
4236c924f48SJed Brown           ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
4246c924f48SJed Brown           if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
425d32f9abdSBarry Smith         }
4267287d2a3SDmitry Karpeev         if (fieldsplit_default || !jac->splitdefined) {
427d32f9abdSBarry Smith           ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
428db4c96c1SJed Brown           for (i=0; i<jac->bs; i++) {
4296c924f48SJed Brown             char splitname[8];
4308caf3d72SBarry Smith             ierr = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr);
4315d4c12cdSJungho Lee             ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr);
43279416396SBarry Smith           }
4335d4c12cdSJungho Lee           jac->defaultsplit = PETSC_TRUE;
434521d7252SBarry Smith         }
43566ffff09SJed Brown       }
4366ce1633cSBarry Smith     }
437edf189efSBarry Smith   } else if (jac->nsplits == 1) {
438edf189efSBarry Smith     if (ilink->is) {
439edf189efSBarry Smith       IS       is2;
440edf189efSBarry Smith       PetscInt nmin,nmax;
441edf189efSBarry Smith 
442edf189efSBarry Smith       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
443edf189efSBarry Smith       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
444db4c96c1SJed Brown       ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr);
445fcfd50ebSBarry Smith       ierr = ISDestroy(&is2);CHKERRQ(ierr);
44682f516ccSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
447edf189efSBarry Smith   }
448d0af7cd3SBarry Smith 
449d0af7cd3SBarry Smith 
45082f516ccSBarry Smith   if (jac->nsplits < 2) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits);
45169a612a9SBarry Smith   PetscFunctionReturn(0);
45269a612a9SBarry Smith }
45369a612a9SBarry Smith 
4545a576424SJed Brown PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(const char pre[], const char name[], char *value[], PetscBool *flg);
455514bf10dSMatthew G Knepley 
45669a612a9SBarry Smith #undef __FUNCT__
45769a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
45869a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
45969a612a9SBarry Smith {
46069a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
46169a612a9SBarry Smith   PetscErrorCode    ierr;
4625a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
4632c9966d7SBarry Smith   PetscInt          i,nsplit;
46469a612a9SBarry Smith   MatStructure      flag = pc->flag;
4655f4ab4e1SJungho Lee   PetscBool         sorted, sorted_col;
46669a612a9SBarry Smith 
46769a612a9SBarry Smith   PetscFunctionBegin;
46869a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
46997bbdb24SBarry Smith   nsplit = jac->nsplits;
4705a9f2f41SSatish Balay   ilink  = jac->head;
47197bbdb24SBarry Smith 
47297bbdb24SBarry Smith   /* get the matrices for each split */
473704ba839SBarry Smith   if (!jac->issetup) {
4741b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
47597bbdb24SBarry Smith 
476704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
477704ba839SBarry Smith 
4785d4c12cdSJungho Lee     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
4792c9966d7SBarry Smith     if (jac->defaultsplit || !ilink->is) {
4802c9966d7SBarry Smith       if (jac->bs <= 0) jac->bs = nsplit;
4812c9966d7SBarry Smith     }
48251f519a2SBarry Smith     bs     = jac->bs;
48397bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
4841b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
4851b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
4861b9fc7fcSBarry Smith       if (jac->defaultsplit) {
487ce94432eSBarry Smith         ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
4885f4ab4e1SJungho Lee         ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr);
489704ba839SBarry Smith       } else if (!ilink->is) {
490ccb205f8SBarry Smith         if (ilink->nfields > 1) {
4915f4ab4e1SJungho Lee           PetscInt *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col;
492785e854fSJed Brown           ierr = PetscMalloc1(ilink->nfields*nslots,&ii);CHKERRQ(ierr);
493785e854fSJed Brown           ierr = PetscMalloc1(ilink->nfields*nslots,&jj);CHKERRQ(ierr);
4941b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
4951b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
4961b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
4975f4ab4e1SJungho Lee               jj[nfields*j + k] = rstart + bs*j + fields_col[k];
49897bbdb24SBarry Smith             }
49997bbdb24SBarry Smith           }
500ce94432eSBarry Smith           ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr);
501ce94432eSBarry Smith           ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr);
502ccb205f8SBarry Smith         } else {
503ce94432eSBarry Smith           ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
504ce94432eSBarry Smith           ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr);
505ccb205f8SBarry Smith        }
5063e197d65SBarry Smith       }
507704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
5085f4ab4e1SJungho Lee       if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); }
5095f4ab4e1SJungho Lee       if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
510704ba839SBarry Smith       ilink = ilink->next;
5111b9fc7fcSBarry Smith     }
5121b9fc7fcSBarry Smith   }
5131b9fc7fcSBarry Smith 
514704ba839SBarry Smith   ilink = jac->head;
51597bbdb24SBarry Smith   if (!jac->pmat) {
516bdddcaaaSMatthew G Knepley     Vec xtmp;
517bdddcaaaSMatthew G Knepley 
5180298fd71SBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,NULL);CHKERRQ(ierr);
519785e854fSJed Brown     ierr = PetscMalloc1(nsplit,&jac->pmat);CHKERRQ(ierr);
520dcca6d9dSJed Brown     ierr = PetscMalloc2(nsplit,&jac->x,nsplit,&jac->y);CHKERRQ(ierr);
521cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
522bdddcaaaSMatthew G Knepley       MatNullSpace sp;
523bdddcaaaSMatthew G Knepley 
524a3df900dSMatthew G Knepley       /* Check for preconditioning matrix attached to IS */
525a3df900dSMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &jac->pmat[i]);CHKERRQ(ierr);
526a3df900dSMatthew G Knepley       if (jac->pmat[i]) {
527a3df900dSMatthew G Knepley         ierr = PetscObjectReference((PetscObject) jac->pmat[i]);CHKERRQ(ierr);
528a3df900dSMatthew G Knepley         if (jac->type == PC_COMPOSITE_SCHUR) {
529a3df900dSMatthew G Knepley           jac->schur_user = jac->pmat[i];
5302fa5cd67SKarl Rupp 
531a3df900dSMatthew G Knepley           ierr = PetscObjectReference((PetscObject) jac->schur_user);CHKERRQ(ierr);
532a3df900dSMatthew G Knepley         }
533a3df900dSMatthew G Knepley       } else {
5345f4ab4e1SJungho Lee         ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
535a3df900dSMatthew G Knepley       }
536bdddcaaaSMatthew G Knepley       /* create work vectors for each split */
537bdddcaaaSMatthew G Knepley       ierr     = MatGetVecs(jac->pmat[i],&jac->x[i],&jac->y[i]);CHKERRQ(ierr);
5380298fd71SBarry Smith       ilink->x = jac->x[i]; ilink->y = jac->y[i]; ilink->z = NULL;
539bdddcaaaSMatthew G Knepley       /* compute scatter contexts needed by multiplicative versions and non-default splits */
5400298fd71SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],NULL,&ilink->sctx);CHKERRQ(ierr);
541ed1f0337SMatthew G Knepley       /* Check for null space attached to IS */
542bafc1b83SMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "nullspace", (PetscObject*) &sp);CHKERRQ(ierr);
543bafc1b83SMatthew G Knepley       if (sp) {
544bafc1b83SMatthew G Knepley         ierr = MatSetNullSpace(jac->pmat[i], sp);CHKERRQ(ierr);
545bafc1b83SMatthew G Knepley       }
546ed1f0337SMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "nearnullspace", (PetscObject*) &sp);CHKERRQ(ierr);
547ed1f0337SMatthew G Knepley       if (sp) {
548ed1f0337SMatthew G Knepley         ierr = MatSetNearNullSpace(jac->pmat[i], sp);CHKERRQ(ierr);
549ed1f0337SMatthew G Knepley       }
550704ba839SBarry Smith       ilink = ilink->next;
551cf502942SBarry Smith     }
552bdddcaaaSMatthew G Knepley     ierr = VecDestroy(&xtmp);CHKERRQ(ierr);
55397bbdb24SBarry Smith   } else {
554cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
555a3df900dSMatthew G Knepley       Mat pmat;
556a3df900dSMatthew G Knepley 
557a3df900dSMatthew G Knepley       /* Check for preconditioning matrix attached to IS */
558a3df900dSMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &pmat);CHKERRQ(ierr);
559a3df900dSMatthew G Knepley       if (!pmat) {
5605f4ab4e1SJungho Lee         ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
561a3df900dSMatthew G Knepley       }
562704ba839SBarry Smith       ilink = ilink->next;
563cf502942SBarry Smith     }
56497bbdb24SBarry Smith   }
565f5236f50SJed Brown   if (pc->useAmat) {
566519d70e2SJed Brown     ilink = jac->head;
567519d70e2SJed Brown     if (!jac->mat) {
568785e854fSJed Brown       ierr = PetscMalloc1(nsplit,&jac->mat);CHKERRQ(ierr);
569519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
5705f4ab4e1SJungho Lee         ierr  = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
571519d70e2SJed Brown         ilink = ilink->next;
572519d70e2SJed Brown       }
573519d70e2SJed Brown     } else {
574519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
5755f4ab4e1SJungho Lee         if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);}
576519d70e2SJed Brown         ilink = ilink->next;
577519d70e2SJed Brown       }
578519d70e2SJed Brown     }
579519d70e2SJed Brown   } else {
580519d70e2SJed Brown     jac->mat = jac->pmat;
581519d70e2SJed Brown   }
58297bbdb24SBarry Smith 
5836c8605c2SJed Brown   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
58468dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
58568dd23aaSBarry Smith     ilink = jac->head;
58668dd23aaSBarry Smith     if (!jac->Afield) {
587785e854fSJed Brown       ierr = PetscMalloc1(nsplit,&jac->Afield);CHKERRQ(ierr);
58868dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
5890298fd71SBarry Smith         ierr  = MatGetSubMatrix(pc->mat,ilink->is,NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
59068dd23aaSBarry Smith         ilink = ilink->next;
59168dd23aaSBarry Smith       }
59268dd23aaSBarry Smith     } else {
59368dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
5940298fd71SBarry Smith         ierr  = MatGetSubMatrix(pc->mat,ilink->is,NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
59568dd23aaSBarry Smith         ilink = ilink->next;
59668dd23aaSBarry Smith       }
59768dd23aaSBarry Smith     }
59868dd23aaSBarry Smith   }
59968dd23aaSBarry Smith 
6003b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
6013b224e63SBarry Smith     IS          ccis;
6024aa3045dSJed Brown     PetscInt    rstart,rend;
603093c86ffSJed Brown     char        lscname[256];
604093c86ffSJed Brown     PetscObject LSC_L;
605ce94432eSBarry Smith 
606ce94432eSBarry Smith     if (nsplit != 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
60768dd23aaSBarry Smith 
608e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
609e6cab6aaSJed Brown     ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
610e6cab6aaSJed Brown 
6113b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
6123b224e63SBarry Smith     if (jac->schur) {
6130298fd71SBarry Smith       KSP kspA = jac->head->ksp, kspInner = NULL, kspUpper = jac->kspupper;
614443836d0SMatthew G Knepley 
615fb3147dbSMatthew G Knepley       ierr  = MatSchurComplementGetKSP(jac->schur, &kspInner);CHKERRQ(ierr);
6163b224e63SBarry Smith       ilink = jac->head;
61749bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
6184aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
619fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
6203b224e63SBarry Smith       ilink = ilink->next;
62149bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
6224aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
623fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
624bee83525SDmitry Karpeev       ierr  = MatSchurComplementUpdateSubMatrices(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],pc->flag);CHKERRQ(ierr);
625a7476a74SDmitry Karpeev       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELFP) {
626a7476a74SDmitry Karpeev 	ierr = MatDestroy(&jac->schurp);CHKERRQ(ierr);
627a7476a74SDmitry Karpeev 	ierr = MatSchurComplementGetPmat(jac->schur,MAT_INITIAL_MATRIX,&jac->schurp);CHKERRQ(ierr);
628a7476a74SDmitry Karpeev       }
629443836d0SMatthew G Knepley       if (kspA != kspInner) {
630443836d0SMatthew G Knepley         ierr = KSPSetOperators(kspA,jac->mat[0],jac->pmat[0],pc->flag);CHKERRQ(ierr);
631443836d0SMatthew G Knepley       }
632443836d0SMatthew G Knepley       if (kspUpper != kspA) {
633443836d0SMatthew G Knepley         ierr = KSPSetOperators(kspUpper,jac->mat[0],jac->pmat[0],pc->flag);CHKERRQ(ierr);
634443836d0SMatthew G Knepley       }
635084e4875SJed Brown       ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr);
6363b224e63SBarry Smith     } else {
637bafc1b83SMatthew G Knepley       const char   *Dprefix;
638470b340bSDmitry Karpeev       char         schurprefix[256], schurmatprefix[256];
639514bf10dSMatthew G Knepley       char         schurtestoption[256];
640bdddcaaaSMatthew G Knepley       MatNullSpace sp;
641514bf10dSMatthew G Knepley       PetscBool    flg;
6423b224e63SBarry Smith 
643a04f6461SBarry Smith       /* extract the A01 and A10 matrices */
6443b224e63SBarry Smith       ilink = jac->head;
64549bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
6464aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
647fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
6483b224e63SBarry Smith       ilink = ilink->next;
64949bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
6504aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
651fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
65220252d06SBarry Smith 
653f5236f50SJed Brown       /* Use mat[0] (diagonal block of Amat) preconditioned by pmat[0] to define Schur complement */
65420252d06SBarry Smith       ierr = MatCreate(((PetscObject)jac->mat[0])->comm,&jac->schur);CHKERRQ(ierr);
65520252d06SBarry Smith       ierr = MatSetType(jac->schur,MATSCHURCOMPLEMENT);CHKERRQ(ierr);
656bee83525SDmitry Karpeev       ierr = MatSchurComplementSetSubMatrices(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);CHKERRQ(ierr);
657470b340bSDmitry Karpeev       ierr = PetscSNPrintf(schurmatprefix, sizeof(schurmatprefix), "%sfieldsplit_%s_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr);
658470b340bSDmitry 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? */
659470b340bSDmitry Karpeev       ierr = MatSetOptionsPrefix(jac->schur,schurmatprefix);CHKERRQ(ierr);
660470b340bSDmitry Karpeev       ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
661bdddcaaaSMatthew G Knepley       ierr = MatGetNullSpace(jac->pmat[1], &sp);CHKERRQ(ierr);
66220252d06SBarry Smith       if (sp) {
66320252d06SBarry Smith         ierr = MatSetNullSpace(jac->schur, sp);CHKERRQ(ierr);
66420252d06SBarry Smith       }
66520252d06SBarry Smith 
66620252d06SBarry Smith       ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_inner_", ilink->splitname);CHKERRQ(ierr);
6670298fd71SBarry Smith       ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);CHKERRQ(ierr);
668514bf10dSMatthew G Knepley       if (flg) {
669514bf10dSMatthew G Knepley         DM  dmInner;
67021635b76SJed Brown         KSP kspInner;
67121635b76SJed Brown 
67221635b76SJed Brown         ierr = MatSchurComplementGetKSP(jac->schur, &kspInner);CHKERRQ(ierr);
67321635b76SJed Brown         ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_inner_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr);
67421635b76SJed Brown         /* Indent this deeper to emphasize the "inner" nature of this solver. */
67521635b76SJed Brown         ierr = PetscObjectIncrementTabLevel((PetscObject)kspInner, (PetscObject) pc, 2);CHKERRQ(ierr);
67621635b76SJed Brown         ierr = KSPSetOptionsPrefix(kspInner, schurprefix);CHKERRQ(ierr);
677514bf10dSMatthew G Knepley 
678514bf10dSMatthew G Knepley         /* Set DM for new solver */
679bafc1b83SMatthew G Knepley         ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr);
68021635b76SJed Brown         ierr = KSPSetDM(kspInner, dmInner);CHKERRQ(ierr);
68121635b76SJed Brown         ierr = KSPSetDMActive(kspInner, PETSC_FALSE);CHKERRQ(ierr);
682514bf10dSMatthew G Knepley       } else {
68321635b76SJed Brown          /* Use the outer solver for the inner solve, but revert the KSPPREONLY from PCFieldSplitSetFields_FieldSplit or
68421635b76SJed Brown           * PCFieldSplitSetIS_FieldSplit. We don't want KSPPREONLY because it makes the Schur complement inexact,
68521635b76SJed Brown           * preventing Schur complement reduction to be an accurate solve. Usually when an iterative solver is used for
68621635b76SJed 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
68721635b76SJed Brown           * GMRES the default. Note that it is also common to use PREONLY for S, in which case S may not be used
68821635b76SJed Brown           * directly, and the user is responsible for setting an inexact method for fieldsplit's A^{-1}. */
68921635b76SJed Brown         ierr = KSPSetType(jac->head->ksp,KSPGMRES);CHKERRQ(ierr);
690514bf10dSMatthew G Knepley         ierr = MatSchurComplementSetKSP(jac->schur,jac->head->ksp);CHKERRQ(ierr);
691bafc1b83SMatthew G Knepley       }
6925a9f2f41SSatish Balay       ierr = KSPSetOperators(jac->head->ksp,jac->mat[0],jac->pmat[0],flag);CHKERRQ(ierr);
6935a9f2f41SSatish Balay       ierr = KSPSetFromOptions(jac->head->ksp);CHKERRQ(ierr);
69497bbdb24SBarry Smith       ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
69597bbdb24SBarry Smith 
696443836d0SMatthew G Knepley       ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_upper_", ilink->splitname);CHKERRQ(ierr);
6970298fd71SBarry Smith       ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);CHKERRQ(ierr);
698443836d0SMatthew G Knepley       if (flg) {
699443836d0SMatthew G Knepley         DM dmInner;
700443836d0SMatthew G Knepley 
701443836d0SMatthew G Knepley         ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_upper_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr);
70282f516ccSBarry Smith         ierr = KSPCreate(PetscObjectComm((PetscObject)pc), &jac->kspupper);CHKERRQ(ierr);
703443836d0SMatthew G Knepley         ierr = KSPSetOptionsPrefix(jac->kspupper, schurprefix);CHKERRQ(ierr);
704443836d0SMatthew G Knepley         ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr);
705443836d0SMatthew G Knepley         ierr = KSPSetDM(jac->kspupper, dmInner);CHKERRQ(ierr);
706443836d0SMatthew G Knepley         ierr = KSPSetDMActive(jac->kspupper, PETSC_FALSE);CHKERRQ(ierr);
707443836d0SMatthew G Knepley         ierr = KSPSetFromOptions(jac->kspupper);CHKERRQ(ierr);
708443836d0SMatthew G Knepley         ierr = KSPSetOperators(jac->kspupper,jac->mat[0],jac->pmat[0],flag);CHKERRQ(ierr);
709443836d0SMatthew G Knepley         ierr = VecDuplicate(jac->head->x, &jac->head->z);CHKERRQ(ierr);
710443836d0SMatthew G Knepley       } else {
711443836d0SMatthew G Knepley         jac->kspupper = jac->head->ksp;
712443836d0SMatthew G Knepley         ierr          = PetscObjectReference((PetscObject) jac->head->ksp);CHKERRQ(ierr);
713443836d0SMatthew G Knepley       }
714443836d0SMatthew G Knepley 
715a7476a74SDmitry Karpeev       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELFP) {
716a7476a74SDmitry Karpeev 	ierr = MatSchurComplementGetPmat(jac->schur,MAT_INITIAL_MATRIX,&jac->schurp);CHKERRQ(ierr);
717a7476a74SDmitry Karpeev       }
718ce94432eSBarry Smith       ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&jac->kspschur);CHKERRQ(ierr);
71997bbdb24SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr);
72020252d06SBarry Smith       ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
72197bbdb24SBarry Smith       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
7227233a360SDmitry Karpeev         PC pcschur;
7237233a360SDmitry Karpeev         ierr = KSPGetPC(jac->kspschur,&pcschur);CHKERRQ(ierr);
7247233a360SDmitry Karpeev         ierr = PCSetType(pcschur,PCNONE);CHKERRQ(ierr);
72597bbdb24SBarry Smith         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
726e74569cdSMatthew G. Knepley       } else if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_FULL) {
727e74569cdSMatthew G. Knepley         ierr = MatSchurComplementComputeExplicitOperator(jac->schur, &jac->schur_user);CHKERRQ(ierr);
72897bbdb24SBarry Smith       }
729e74569cdSMatthew G. Knepley       ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
73097bbdb24SBarry Smith       ierr = KSPGetOptionsPrefix(jac->head->next->ksp, &Dprefix);CHKERRQ(ierr);
73197bbdb24SBarry Smith       ierr = KSPSetOptionsPrefix(jac->kspschur,         Dprefix);CHKERRQ(ierr);
732b20b4189SMatthew G. Knepley       /* propogate DM */
733b20b4189SMatthew G. Knepley       {
734b20b4189SMatthew G. Knepley         DM sdm;
735b20b4189SMatthew G. Knepley         ierr = KSPGetDM(jac->head->next->ksp, &sdm);CHKERRQ(ierr);
736b20b4189SMatthew G. Knepley         if (sdm) {
737b20b4189SMatthew G. Knepley           ierr = KSPSetDM(jac->kspschur, sdm);CHKERRQ(ierr);
738b20b4189SMatthew G. Knepley           ierr = KSPSetDMActive(jac->kspschur, PETSC_FALSE);CHKERRQ(ierr);
739b20b4189SMatthew G. Knepley         }
740b20b4189SMatthew G. Knepley       }
74197bbdb24SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
74297bbdb24SBarry Smith       /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */
74397bbdb24SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
74497bbdb24SBarry Smith     }
74597bbdb24SBarry Smith 
7465a9f2f41SSatish Balay     /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */
7478caf3d72SBarry Smith     ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_L",ilink->splitname);CHKERRQ(ierr);
748519d70e2SJed Brown     ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
7493b224e63SBarry Smith     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
750c1570756SJed Brown     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);}
7518caf3d72SBarry Smith     ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr);
75297bbdb24SBarry Smith     ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
7535a9f2f41SSatish Balay     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
7540971522cSBarry Smith     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);}
75597bbdb24SBarry Smith   } else {
75668bd789dSDmitry Karpeev     /* set up the individual splits' PCs */
7570971522cSBarry Smith     i     = 0;
7580971522cSBarry Smith     ilink = jac->head;
7590971522cSBarry Smith     while (ilink) {
7600971522cSBarry Smith       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr);
7610971522cSBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
7620971522cSBarry Smith       if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);}
7630971522cSBarry Smith       i++;
7640971522cSBarry Smith       ilink = ilink->next;
7650971522cSBarry Smith     }
7663b224e63SBarry Smith   }
7673b224e63SBarry Smith 
768c1570756SJed Brown   jac->suboptionsset = PETSC_TRUE;
7690971522cSBarry Smith   PetscFunctionReturn(0);
7700971522cSBarry Smith }
7710971522cSBarry Smith 
7725a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
773ca9f406cSSatish Balay   (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
774ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
7755a9f2f41SSatish Balay    KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
776ca9f406cSSatish Balay    VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
777ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
77879416396SBarry Smith 
7790971522cSBarry Smith #undef __FUNCT__
7803b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur"
7813b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
7823b224e63SBarry Smith {
7833b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
7843b224e63SBarry Smith   PetscErrorCode    ierr;
7853b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
786443836d0SMatthew G Knepley   KSP               kspA   = ilinkA->ksp, kspLower = kspA, kspUpper = jac->kspupper;
7873b224e63SBarry Smith 
7883b224e63SBarry Smith   PetscFunctionBegin;
789c5d2311dSJed Brown   switch (jac->schurfactorization) {
790c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_DIAG:
791a04f6461SBarry Smith     /* [A00 0; 0 -S], positive definite, suitable for MINRES */
792c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
793c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
794c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
795443836d0SMatthew G Knepley     ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
796c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
797c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
798c5d2311dSJed Brown     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
799c5d2311dSJed Brown     ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr);
800c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
801c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
802c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
803c5d2311dSJed Brown     break;
804c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_LOWER:
805a04f6461SBarry Smith     /* [A00 0; A10 S], suitable for left preconditioning */
806c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
807c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
808443836d0SMatthew G Knepley     ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
809c5d2311dSJed Brown     ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
810c5d2311dSJed Brown     ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr);
811c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
812c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
813c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
814c5d2311dSJed Brown     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
815c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
816c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
817c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
818c5d2311dSJed Brown     break;
819c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_UPPER:
820a04f6461SBarry Smith     /* [A00 A01; 0 S], suitable for right preconditioning */
821c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
822c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
823c5d2311dSJed Brown     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
824c5d2311dSJed Brown     ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
825c5d2311dSJed Brown     ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr);
826c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
827c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
828c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
829443836d0SMatthew G Knepley     ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
830c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
831c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
832c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
833c5d2311dSJed Brown     break;
834c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_FULL:
835a04f6461SBarry 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 */
8363b224e63SBarry Smith     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8373b224e63SBarry Smith     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
838443836d0SMatthew G Knepley     ierr = KSPSolve(kspLower,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
8393b224e63SBarry Smith     ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
8403b224e63SBarry Smith     ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
8413b224e63SBarry Smith     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8423b224e63SBarry Smith     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8433b224e63SBarry Smith 
8443b224e63SBarry Smith     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
8453b224e63SBarry Smith     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8463b224e63SBarry Smith     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8473b224e63SBarry Smith 
848443836d0SMatthew G Knepley     if (kspUpper == kspA) {
8493b224e63SBarry Smith       ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
8503b224e63SBarry Smith       ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
851443836d0SMatthew G Knepley       ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
852443836d0SMatthew G Knepley     } else {
853443836d0SMatthew G Knepley       ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
854443836d0SMatthew G Knepley       ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
855443836d0SMatthew G Knepley       ierr = KSPSolve(kspUpper,ilinkA->x,ilinkA->z);CHKERRQ(ierr);
856443836d0SMatthew G Knepley       ierr = VecAXPY(ilinkA->y,-1.0,ilinkA->z);CHKERRQ(ierr);
857443836d0SMatthew G Knepley     }
8583b224e63SBarry Smith     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8593b224e63SBarry Smith     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
860c5d2311dSJed Brown   }
8613b224e63SBarry Smith   PetscFunctionReturn(0);
8623b224e63SBarry Smith }
8633b224e63SBarry Smith 
8643b224e63SBarry Smith #undef __FUNCT__
8650971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
8660971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
8670971522cSBarry Smith {
8680971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
8690971522cSBarry Smith   PetscErrorCode    ierr;
8705a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
871939b8a20SBarry Smith   PetscInt          cnt,bs;
8720971522cSBarry Smith 
8730971522cSBarry Smith   PetscFunctionBegin;
87479416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
8751b9fc7fcSBarry Smith     if (jac->defaultsplit) {
876939b8a20SBarry Smith       ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
877ce94432eSBarry 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);
878939b8a20SBarry Smith       ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr);
879ce94432eSBarry 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);
8800971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
8815a9f2f41SSatish Balay       while (ilink) {
8825a9f2f41SSatish Balay         ierr  = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
8835a9f2f41SSatish Balay         ilink = ilink->next;
8840971522cSBarry Smith       }
8850971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
8861b9fc7fcSBarry Smith     } else {
887efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
8885a9f2f41SSatish Balay       while (ilink) {
8895a9f2f41SSatish Balay         ierr  = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
8905a9f2f41SSatish Balay         ilink = ilink->next;
8911b9fc7fcSBarry Smith       }
8921b9fc7fcSBarry Smith     }
89316913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
89479416396SBarry Smith     if (!jac->w1) {
89579416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
89679416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
89779416396SBarry Smith     }
898efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
8995a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
9003e197d65SBarry Smith     cnt  = 1;
9015a9f2f41SSatish Balay     while (ilink->next) {
9025a9f2f41SSatish Balay       ilink = ilink->next;
9033e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
9043e197d65SBarry Smith       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
9053e197d65SBarry Smith       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
9063e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9073e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9083e197d65SBarry Smith       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
9093e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9103e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9113e197d65SBarry Smith     }
91251f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
91311755939SBarry Smith       cnt -= 2;
91451f519a2SBarry Smith       while (ilink->previous) {
91551f519a2SBarry Smith         ilink = ilink->previous;
91611755939SBarry Smith         /* compute the residual only over the part of the vector needed */
91711755939SBarry Smith         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
91811755939SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
91911755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
92011755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
92111755939SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
92211755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
92311755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
92451f519a2SBarry Smith       }
92511755939SBarry Smith     }
926ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
9270971522cSBarry Smith   PetscFunctionReturn(0);
9280971522cSBarry Smith }
9290971522cSBarry Smith 
930421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
931ca9f406cSSatish Balay   (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
932ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
933421e10b8SBarry Smith    KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
934ca9f406cSSatish Balay    VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
935ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
936421e10b8SBarry Smith 
937421e10b8SBarry Smith #undef __FUNCT__
9388c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit"
939421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
940421e10b8SBarry Smith {
941421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
942421e10b8SBarry Smith   PetscErrorCode    ierr;
943421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
944939b8a20SBarry Smith   PetscInt          bs;
945421e10b8SBarry Smith 
946421e10b8SBarry Smith   PetscFunctionBegin;
947421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
948421e10b8SBarry Smith     if (jac->defaultsplit) {
949939b8a20SBarry Smith       ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
950ce94432eSBarry 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);
951939b8a20SBarry Smith       ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr);
952ce94432eSBarry 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);
953421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
954421e10b8SBarry Smith       while (ilink) {
955421e10b8SBarry Smith         ierr  = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
956421e10b8SBarry Smith         ilink = ilink->next;
957421e10b8SBarry Smith       }
958421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
959421e10b8SBarry Smith     } else {
960421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
961421e10b8SBarry Smith       while (ilink) {
962421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
963421e10b8SBarry Smith         ilink = ilink->next;
964421e10b8SBarry Smith       }
965421e10b8SBarry Smith     }
966421e10b8SBarry Smith   } else {
967421e10b8SBarry Smith     if (!jac->w1) {
968421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
969421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
970421e10b8SBarry Smith     }
971421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
972421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
973421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
974421e10b8SBarry Smith       while (ilink->next) {
975421e10b8SBarry Smith         ilink = ilink->next;
9769989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
977421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
978421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
979421e10b8SBarry Smith       }
980421e10b8SBarry Smith       while (ilink->previous) {
981421e10b8SBarry Smith         ilink = ilink->previous;
9829989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
983421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
984421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
985421e10b8SBarry Smith       }
986421e10b8SBarry Smith     } else {
987421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
988421e10b8SBarry Smith         ilink = ilink->next;
989421e10b8SBarry Smith       }
990421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
991421e10b8SBarry Smith       while (ilink->previous) {
992421e10b8SBarry Smith         ilink = ilink->previous;
9939989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
994421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
995421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
996421e10b8SBarry Smith       }
997421e10b8SBarry Smith     }
998421e10b8SBarry Smith   }
999421e10b8SBarry Smith   PetscFunctionReturn(0);
1000421e10b8SBarry Smith }
1001421e10b8SBarry Smith 
10020971522cSBarry Smith #undef __FUNCT__
1003574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit"
1004574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc)
10050971522cSBarry Smith {
10060971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
10070971522cSBarry Smith   PetscErrorCode    ierr;
10085a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
10090971522cSBarry Smith 
10100971522cSBarry Smith   PetscFunctionBegin;
10115a9f2f41SSatish Balay   while (ilink) {
1012574deadeSBarry Smith     ierr  = KSPReset(ilink->ksp);CHKERRQ(ierr);
1013fcfd50ebSBarry Smith     ierr  = VecDestroy(&ilink->x);CHKERRQ(ierr);
1014fcfd50ebSBarry Smith     ierr  = VecDestroy(&ilink->y);CHKERRQ(ierr);
1015443836d0SMatthew G Knepley     ierr  = VecDestroy(&ilink->z);CHKERRQ(ierr);
1016fcfd50ebSBarry Smith     ierr  = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr);
1017fcfd50ebSBarry Smith     ierr  = ISDestroy(&ilink->is);CHKERRQ(ierr);
1018b5787286SJed Brown     ierr  = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
10195a9f2f41SSatish Balay     next  = ilink->next;
10205a9f2f41SSatish Balay     ilink = next;
10210971522cSBarry Smith   }
102205b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
1023574deadeSBarry Smith   if (jac->mat && jac->mat != jac->pmat) {
1024574deadeSBarry Smith     ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);
1025574deadeSBarry Smith   } else if (jac->mat) {
10260298fd71SBarry Smith     jac->mat = NULL;
1027574deadeSBarry Smith   }
102897bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
102968dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
10306bf464f9SBarry Smith   ierr       = VecDestroy(&jac->w1);CHKERRQ(ierr);
10316bf464f9SBarry Smith   ierr       = VecDestroy(&jac->w2);CHKERRQ(ierr);
10326bf464f9SBarry Smith   ierr       = MatDestroy(&jac->schur);CHKERRQ(ierr);
1033a7476a74SDmitry Karpeev   ierr       = MatDestroy(&jac->schurp);CHKERRQ(ierr);
10346bf464f9SBarry Smith   ierr       = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
10356bf464f9SBarry Smith   ierr       = KSPDestroy(&jac->kspschur);CHKERRQ(ierr);
1036d78dad28SBarry Smith   ierr       = KSPDestroy(&jac->kspupper);CHKERRQ(ierr);
10376bf464f9SBarry Smith   ierr       = MatDestroy(&jac->B);CHKERRQ(ierr);
10386bf464f9SBarry Smith   ierr       = MatDestroy(&jac->C);CHKERRQ(ierr);
103963ec74ffSBarry Smith   jac->reset = PETSC_TRUE;
1040574deadeSBarry Smith   PetscFunctionReturn(0);
1041574deadeSBarry Smith }
1042574deadeSBarry Smith 
1043574deadeSBarry Smith #undef __FUNCT__
1044574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
1045574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
1046574deadeSBarry Smith {
1047574deadeSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1048574deadeSBarry Smith   PetscErrorCode    ierr;
1049574deadeSBarry Smith   PC_FieldSplitLink ilink = jac->head,next;
1050574deadeSBarry Smith 
1051574deadeSBarry Smith   PetscFunctionBegin;
1052574deadeSBarry Smith   ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr);
1053574deadeSBarry Smith   while (ilink) {
10546bf464f9SBarry Smith     ierr  = KSPDestroy(&ilink->ksp);CHKERRQ(ierr);
1055574deadeSBarry Smith     next  = ilink->next;
1056574deadeSBarry Smith     ierr  = PetscFree(ilink->splitname);CHKERRQ(ierr);
1057574deadeSBarry Smith     ierr  = PetscFree(ilink->fields);CHKERRQ(ierr);
10585d4c12cdSJungho Lee     ierr  = PetscFree(ilink->fields_col);CHKERRQ(ierr);
1059574deadeSBarry Smith     ierr  = PetscFree(ilink);CHKERRQ(ierr);
1060574deadeSBarry Smith     ilink = next;
1061574deadeSBarry Smith   }
1062574deadeSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
1063c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1064bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",NULL);CHKERRQ(ierr);
1065bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",NULL);CHKERRQ(ierr);
1066bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",NULL);CHKERRQ(ierr);
1067bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",NULL);CHKERRQ(ierr);
1068bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",NULL);CHKERRQ(ierr);
1069bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",NULL);CHKERRQ(ierr);
1070bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",NULL);CHKERRQ(ierr);
10710971522cSBarry Smith   PetscFunctionReturn(0);
10720971522cSBarry Smith }
10730971522cSBarry Smith 
10740971522cSBarry Smith #undef __FUNCT__
10750971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
10760971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
10770971522cSBarry Smith {
10781b9fc7fcSBarry Smith   PetscErrorCode  ierr;
10796c924f48SJed Brown   PetscInt        bs;
1080bc59fbc5SBarry Smith   PetscBool       flg,stokes = PETSC_FALSE;
10819dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
10823b224e63SBarry Smith   PCCompositeType ctype;
10831b9fc7fcSBarry Smith 
10840971522cSBarry Smith   PetscFunctionBegin;
10851b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
10864ab8060aSDmitry Karpeev   ierr = PetscOptionsBool("-pc_fieldsplit_dm_splits","Whether to use DMCreateFieldDecomposition() for splits","PCFieldSplitSetDMSplits",jac->dm_splits,&jac->dm_splits,NULL);CHKERRQ(ierr);
108751f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
108851f519a2SBarry Smith   if (flg) {
108951f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
109051f519a2SBarry Smith   }
1091704ba839SBarry Smith 
10920298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);CHKERRQ(ierr);
1093c0adfefeSBarry Smith   if (stokes) {
1094c0adfefeSBarry Smith     ierr          = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
1095c0adfefeSBarry Smith     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
1096c0adfefeSBarry Smith   }
1097c0adfefeSBarry Smith 
10983b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
10993b224e63SBarry Smith   if (flg) {
11003b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
11013b224e63SBarry Smith   }
1102c30613efSMatthew Knepley   /* Only setup fields once */
1103c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
1104d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
1105d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
11066c924f48SJed Brown     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
11076c924f48SJed Brown     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
1108d32f9abdSBarry Smith   }
1109c5d2311dSJed Brown   if (jac->type == PC_COMPOSITE_SCHUR) {
1110c9c6ffaaSJed Brown     ierr = PetscOptionsGetEnum(((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr);
1111c9c6ffaaSJed Brown     if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);}
11120298fd71SBarry 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);
11130298fd71SBarry Smith     ierr = PetscOptionsEnum("-pc_fieldsplit_schur_precondition","How to build preconditioner for Schur complement","PCFieldSplitSchurPrecondition",PCFieldSplitSchurPreTypes,(PetscEnum)jac->schurpre,(PetscEnum*)&jac->schurpre,NULL);CHKERRQ(ierr);
1114c5d2311dSJed Brown   }
11151b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
11160971522cSBarry Smith   PetscFunctionReturn(0);
11170971522cSBarry Smith }
11180971522cSBarry Smith 
11190971522cSBarry Smith /*------------------------------------------------------------------------------------*/
11200971522cSBarry Smith 
11210971522cSBarry Smith #undef __FUNCT__
11220971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
11231e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
11240971522cSBarry Smith {
112597bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
11260971522cSBarry Smith   PetscErrorCode    ierr;
11275a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
112869a612a9SBarry Smith   char              prefix[128];
11295d4c12cdSJungho Lee   PetscInt          i;
11300971522cSBarry Smith 
11310971522cSBarry Smith   PetscFunctionBegin;
11326c924f48SJed Brown   if (jac->splitdefined) {
11336c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
11346c924f48SJed Brown     PetscFunctionReturn(0);
11356c924f48SJed Brown   }
113651f519a2SBarry Smith   for (i=0; i<n; i++) {
1137e32f2f54SBarry 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);
1138e32f2f54SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
113951f519a2SBarry Smith   }
1140b00a9115SJed Brown   ierr = PetscNew(&ilink);CHKERRQ(ierr);
1141a04f6461SBarry Smith   if (splitname) {
1142db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1143a04f6461SBarry Smith   } else {
1144785e854fSJed Brown     ierr = PetscMalloc1(3,&ilink->splitname);CHKERRQ(ierr);
1145a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
1146a04f6461SBarry Smith   }
1147785e854fSJed Brown   ierr = PetscMalloc1(n,&ilink->fields);CHKERRQ(ierr);
11485a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
1149785e854fSJed Brown   ierr = PetscMalloc1(n,&ilink->fields_col);CHKERRQ(ierr);
11505d4c12cdSJungho Lee   ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr);
11512fa5cd67SKarl Rupp 
11525a9f2f41SSatish Balay   ilink->nfields = n;
11530298fd71SBarry Smith   ilink->next    = NULL;
1154ce94432eSBarry Smith   ierr           = KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);CHKERRQ(ierr);
115520252d06SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
11565a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
11579005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
115869a612a9SBarry Smith 
11598caf3d72SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);CHKERRQ(ierr);
11605a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
11610971522cSBarry Smith 
11620971522cSBarry Smith   if (!next) {
11635a9f2f41SSatish Balay     jac->head       = ilink;
11640298fd71SBarry Smith     ilink->previous = NULL;
11650971522cSBarry Smith   } else {
11660971522cSBarry Smith     while (next->next) {
11670971522cSBarry Smith       next = next->next;
11680971522cSBarry Smith     }
11695a9f2f41SSatish Balay     next->next      = ilink;
117051f519a2SBarry Smith     ilink->previous = next;
11710971522cSBarry Smith   }
11720971522cSBarry Smith   jac->nsplits++;
11730971522cSBarry Smith   PetscFunctionReturn(0);
11740971522cSBarry Smith }
11750971522cSBarry Smith 
1176e69d4d44SBarry Smith #undef __FUNCT__
1177e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
11781e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
1179e69d4d44SBarry Smith {
1180e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1181e69d4d44SBarry Smith   PetscErrorCode ierr;
1182e69d4d44SBarry Smith 
1183e69d4d44SBarry Smith   PetscFunctionBegin;
1184785e854fSJed Brown   ierr = PetscMalloc1(jac->nsplits,subksp);CHKERRQ(ierr);
1185e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
11862fa5cd67SKarl Rupp 
1187e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
118813e0d083SBarry Smith   if (n) *n = jac->nsplits;
1189e69d4d44SBarry Smith   PetscFunctionReturn(0);
1190e69d4d44SBarry Smith }
11910971522cSBarry Smith 
11920971522cSBarry Smith #undef __FUNCT__
119369a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
11941e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
11950971522cSBarry Smith {
11960971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
11970971522cSBarry Smith   PetscErrorCode    ierr;
11980971522cSBarry Smith   PetscInt          cnt   = 0;
11995a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
12000971522cSBarry Smith 
12010971522cSBarry Smith   PetscFunctionBegin;
1202785e854fSJed Brown   ierr = PetscMalloc1(jac->nsplits,subksp);CHKERRQ(ierr);
12035a9f2f41SSatish Balay   while (ilink) {
12045a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
12055a9f2f41SSatish Balay     ilink            = ilink->next;
12060971522cSBarry Smith   }
12075d480477SMatthew 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);
120813e0d083SBarry Smith   if (n) *n = jac->nsplits;
12090971522cSBarry Smith   PetscFunctionReturn(0);
12100971522cSBarry Smith }
12110971522cSBarry Smith 
1212704ba839SBarry Smith #undef __FUNCT__
1213704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
12141e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
1215704ba839SBarry Smith {
1216704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1217704ba839SBarry Smith   PetscErrorCode    ierr;
1218704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
1219704ba839SBarry Smith   char              prefix[128];
1220704ba839SBarry Smith 
1221704ba839SBarry Smith   PetscFunctionBegin;
12226c924f48SJed Brown   if (jac->splitdefined) {
12236c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
12246c924f48SJed Brown     PetscFunctionReturn(0);
12256c924f48SJed Brown   }
1226b00a9115SJed Brown   ierr = PetscNew(&ilink);CHKERRQ(ierr);
1227a04f6461SBarry Smith   if (splitname) {
1228db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1229a04f6461SBarry Smith   } else {
1230785e854fSJed Brown     ierr = PetscMalloc1(8,&ilink->splitname);CHKERRQ(ierr);
1231b5787286SJed Brown     ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr);
1232a04f6461SBarry Smith   }
12331ab39975SBarry Smith   ierr          = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1234b5787286SJed Brown   ierr          = ISDestroy(&ilink->is);CHKERRQ(ierr);
1235b5787286SJed Brown   ilink->is     = is;
1236b5787286SJed Brown   ierr          = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1237b5787286SJed Brown   ierr          = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
1238b5787286SJed Brown   ilink->is_col = is;
12390298fd71SBarry Smith   ilink->next   = NULL;
1240ce94432eSBarry Smith   ierr          = KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);CHKERRQ(ierr);
124120252d06SBarry Smith   ierr          = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1242704ba839SBarry Smith   ierr          = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
12439005cf84SBarry Smith   ierr          = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
1244704ba839SBarry Smith 
12458caf3d72SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);CHKERRQ(ierr);
1246704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
1247704ba839SBarry Smith 
1248704ba839SBarry Smith   if (!next) {
1249704ba839SBarry Smith     jac->head       = ilink;
12500298fd71SBarry Smith     ilink->previous = NULL;
1251704ba839SBarry Smith   } else {
1252704ba839SBarry Smith     while (next->next) {
1253704ba839SBarry Smith       next = next->next;
1254704ba839SBarry Smith     }
1255704ba839SBarry Smith     next->next      = ilink;
1256704ba839SBarry Smith     ilink->previous = next;
1257704ba839SBarry Smith   }
1258704ba839SBarry Smith   jac->nsplits++;
1259704ba839SBarry Smith   PetscFunctionReturn(0);
1260704ba839SBarry Smith }
1261704ba839SBarry Smith 
12620971522cSBarry Smith #undef __FUNCT__
12630971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
12640971522cSBarry Smith /*@
12650971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
12660971522cSBarry Smith 
1267ad4df100SBarry Smith     Logically Collective on PC
12680971522cSBarry Smith 
12690971522cSBarry Smith     Input Parameters:
12700971522cSBarry Smith +   pc  - the preconditioner context
12710298fd71SBarry Smith .   splitname - name of this split, if NULL the number of the split is used
12720971522cSBarry Smith .   n - the number of fields in this split
1273db4c96c1SJed Brown -   fields - the fields in this split
12740971522cSBarry Smith 
12750971522cSBarry Smith     Level: intermediate
12760971522cSBarry Smith 
1277d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
1278d32f9abdSBarry Smith 
12797287d2a3SDmitry Karpeev      The PCFieldSplitSetFields() is for defining fields as strided blocks. For example, if the block
1280d32f9abdSBarry 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
1281d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1282d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
1283d32f9abdSBarry Smith 
1284db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
1285db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1286db4c96c1SJed Brown 
12875d4c12cdSJungho Lee      Developer Note: This routine does not actually create the IS representing the split, that is delayed
12885d4c12cdSJungho Lee      until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be
12895d4c12cdSJungho Lee      available when this routine is called.
12905d4c12cdSJungho Lee 
1291d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
12920971522cSBarry Smith 
12930971522cSBarry Smith @*/
12945d4c12cdSJungho Lee PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
12950971522cSBarry Smith {
12964ac538c5SBarry Smith   PetscErrorCode ierr;
12970971522cSBarry Smith 
12980971522cSBarry Smith   PetscFunctionBegin;
12990700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1300db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1301ce94432eSBarry Smith   if (n < 1) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1302db4c96c1SJed Brown   PetscValidIntPointer(fields,3);
13035d4c12cdSJungho Lee   ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt*,const PetscInt*),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr);
13040971522cSBarry Smith   PetscFunctionReturn(0);
13050971522cSBarry Smith }
13060971522cSBarry Smith 
13070971522cSBarry Smith #undef __FUNCT__
1308704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
1309704ba839SBarry Smith /*@
1310704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
1311704ba839SBarry Smith 
1312ad4df100SBarry Smith     Logically Collective on PC
1313704ba839SBarry Smith 
1314704ba839SBarry Smith     Input Parameters:
1315704ba839SBarry Smith +   pc  - the preconditioner context
13160298fd71SBarry Smith .   splitname - name of this split, if NULL the number of the split is used
1317db4c96c1SJed Brown -   is - the index set that defines the vector elements in this field
1318704ba839SBarry Smith 
1319d32f9abdSBarry Smith 
1320a6ffb8dbSJed Brown     Notes:
1321a6ffb8dbSJed Brown     Use PCFieldSplitSetFields(), for fields defined by strided types.
1322a6ffb8dbSJed Brown 
1323db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
1324db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1325d32f9abdSBarry Smith 
1326704ba839SBarry Smith     Level: intermediate
1327704ba839SBarry Smith 
1328704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
1329704ba839SBarry Smith 
1330704ba839SBarry Smith @*/
13317087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1332704ba839SBarry Smith {
13334ac538c5SBarry Smith   PetscErrorCode ierr;
1334704ba839SBarry Smith 
1335704ba839SBarry Smith   PetscFunctionBegin;
13360700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13377b62db95SJungho Lee   if (splitname) PetscValidCharPointer(splitname,2);
1338db4c96c1SJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,3);
13394ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr);
1340704ba839SBarry Smith   PetscFunctionReturn(0);
1341704ba839SBarry Smith }
1342704ba839SBarry Smith 
1343704ba839SBarry Smith #undef __FUNCT__
134457a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS"
134557a9adfeSMatthew G Knepley /*@
134657a9adfeSMatthew G Knepley     PCFieldSplitGetIS - Retrieves the elements for a field as an IS
134757a9adfeSMatthew G Knepley 
134857a9adfeSMatthew G Knepley     Logically Collective on PC
134957a9adfeSMatthew G Knepley 
135057a9adfeSMatthew G Knepley     Input Parameters:
135157a9adfeSMatthew G Knepley +   pc  - the preconditioner context
135257a9adfeSMatthew G Knepley -   splitname - name of this split
135357a9adfeSMatthew G Knepley 
135457a9adfeSMatthew G Knepley     Output Parameter:
13550298fd71SBarry Smith -   is - the index set that defines the vector elements in this field, or NULL if the field is not found
135657a9adfeSMatthew G Knepley 
135757a9adfeSMatthew G Knepley     Level: intermediate
135857a9adfeSMatthew G Knepley 
135957a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS()
136057a9adfeSMatthew G Knepley 
136157a9adfeSMatthew G Knepley @*/
136257a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
136357a9adfeSMatthew G Knepley {
136457a9adfeSMatthew G Knepley   PetscErrorCode ierr;
136557a9adfeSMatthew G Knepley 
136657a9adfeSMatthew G Knepley   PetscFunctionBegin;
136757a9adfeSMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
136857a9adfeSMatthew G Knepley   PetscValidCharPointer(splitname,2);
136957a9adfeSMatthew G Knepley   PetscValidPointer(is,3);
137057a9adfeSMatthew G Knepley   {
137157a9adfeSMatthew G Knepley     PC_FieldSplit     *jac  = (PC_FieldSplit*) pc->data;
137257a9adfeSMatthew G Knepley     PC_FieldSplitLink ilink = jac->head;
137357a9adfeSMatthew G Knepley     PetscBool         found;
137457a9adfeSMatthew G Knepley 
13750298fd71SBarry Smith     *is = NULL;
137657a9adfeSMatthew G Knepley     while (ilink) {
137757a9adfeSMatthew G Knepley       ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr);
137857a9adfeSMatthew G Knepley       if (found) {
137957a9adfeSMatthew G Knepley         *is = ilink->is;
138057a9adfeSMatthew G Knepley         break;
138157a9adfeSMatthew G Knepley       }
138257a9adfeSMatthew G Knepley       ilink = ilink->next;
138357a9adfeSMatthew G Knepley     }
138457a9adfeSMatthew G Knepley   }
138557a9adfeSMatthew G Knepley   PetscFunctionReturn(0);
138657a9adfeSMatthew G Knepley }
138757a9adfeSMatthew G Knepley 
138857a9adfeSMatthew G Knepley #undef __FUNCT__
138951f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
139051f519a2SBarry Smith /*@
139151f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
139251f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
139351f519a2SBarry Smith 
1394ad4df100SBarry Smith     Logically Collective on PC
139551f519a2SBarry Smith 
139651f519a2SBarry Smith     Input Parameters:
139751f519a2SBarry Smith +   pc  - the preconditioner context
139851f519a2SBarry Smith -   bs - the block size
139951f519a2SBarry Smith 
140051f519a2SBarry Smith     Level: intermediate
140151f519a2SBarry Smith 
140251f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
140351f519a2SBarry Smith 
140451f519a2SBarry Smith @*/
14057087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
140651f519a2SBarry Smith {
14074ac538c5SBarry Smith   PetscErrorCode ierr;
140851f519a2SBarry Smith 
140951f519a2SBarry Smith   PetscFunctionBegin;
14100700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1411c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,bs,2);
14124ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr);
141351f519a2SBarry Smith   PetscFunctionReturn(0);
141451f519a2SBarry Smith }
141551f519a2SBarry Smith 
141651f519a2SBarry Smith #undef __FUNCT__
141769a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
14180971522cSBarry Smith /*@C
141969a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
14200971522cSBarry Smith 
142169a612a9SBarry Smith    Collective on KSP
14220971522cSBarry Smith 
14230971522cSBarry Smith    Input Parameter:
14240971522cSBarry Smith .  pc - the preconditioner context
14250971522cSBarry Smith 
14260971522cSBarry Smith    Output Parameters:
142713e0d083SBarry Smith +  n - the number of splits
142869a612a9SBarry Smith -  pc - the array of KSP contexts
14290971522cSBarry Smith 
14300971522cSBarry Smith    Note:
1431d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
1432d32f9abdSBarry Smith    (not the KSP just the array that contains them).
14330971522cSBarry Smith 
143469a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
14350971522cSBarry Smith 
1436196cc216SBarry Smith    Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs.
14370298fd71SBarry Smith       You can call PCFieldSplitGetSubKSP(pc,n,NULL_OBJECT,ierr) to determine how large the
1438196cc216SBarry Smith       KSP array must be.
1439196cc216SBarry Smith 
1440196cc216SBarry Smith 
14410971522cSBarry Smith    Level: advanced
14420971522cSBarry Smith 
14430971522cSBarry Smith .seealso: PCFIELDSPLIT
14440971522cSBarry Smith @*/
14457087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
14460971522cSBarry Smith {
14474ac538c5SBarry Smith   PetscErrorCode ierr;
14480971522cSBarry Smith 
14490971522cSBarry Smith   PetscFunctionBegin;
14500700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
145113e0d083SBarry Smith   if (n) PetscValidIntPointer(n,2);
14524ac538c5SBarry Smith   ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
14530971522cSBarry Smith   PetscFunctionReturn(0);
14540971522cSBarry Smith }
14550971522cSBarry Smith 
1456e69d4d44SBarry Smith #undef __FUNCT__
1457e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition"
1458e69d4d44SBarry Smith /*@
1459e69d4d44SBarry Smith     PCFieldSplitSchurPrecondition -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
1460a04f6461SBarry Smith       A11 matrix. Otherwise no preconditioner is used.
1461e69d4d44SBarry Smith 
1462e69d4d44SBarry Smith     Collective on PC
1463e69d4d44SBarry Smith 
1464e69d4d44SBarry Smith     Input Parameters:
1465e69d4d44SBarry Smith +   pc      - the preconditioner context
1466a7476a74SDmitry Karpeev .   ptype   - which matrix to use for preconditioning the Schur complement: PC_FIELDSPLIT_SCHUR_PRE_A11 (default), PC_FIELDSPLIT_SCHUR_PRE_SELF, PC_FIELDSPLIT_PRE_USER
14670298fd71SBarry Smith -   userpre - matrix to use for preconditioning, or NULL
1468084e4875SJed Brown 
1469e69d4d44SBarry Smith     Options Database:
1470*2e71c61dSMatthew G. Knepley .     -pc_fieldsplit_schur_precondition <self,selfp,user,a11,full> default is a11
1471e69d4d44SBarry Smith 
1472fd1303e9SJungho Lee     Notes:
1473fd1303e9SJungho Lee $    If ptype is
1474fd1303e9SJungho Lee $        user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument
1475fd1303e9SJungho Lee $             to this function).
1476e87fab1bSBarry Smith $        a11 then the preconditioner for the Schur complement is generated by the block diagonal part of the original
1477fd1303e9SJungho Lee $             matrix associated with the Schur complement (i.e. A11)
1478e74569cdSMatthew G. Knepley $        full then the preconditioner uses the exact Schur complement (this is expensive)
1479fd1303e9SJungho Lee $        self the preconditioner for the Schur complement is generated from the Schur complement matrix itself:
1480fd1303e9SJungho Lee $             The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC
1481fd1303e9SJungho Lee $             preconditioner
1482a7476a74SDmitry Karpeev $        selfp then the preconditioning matrix is an explicitly-assembled approximation Sp = A11 - A10 inv(diag(A00)) A01
1483a7476a74SDmitry Karpeev $             This is only a good preconditioner when diag(A00) is a good preconditioner for A00.
1484fd1303e9SJungho Lee 
1485e87fab1bSBarry Smith      When solving a saddle point problem, where the A11 block is identically zero, using a11 as the ptype only makes sense
1486fd1303e9SJungho Lee     with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
1487a7476a74SDmitry Karpeev     -fieldsplit_1_pc_type lsc which uses the least squares commutator to compute a preconditioner for the Schur complement.
1488fd1303e9SJungho Lee 
1489a7476a74SDmitry Karpeev     Developer Notes: This is a terrible name, which gives no good indication of what the function does.
1490a7476a74SDmitry Karpeev                      Among other things, it should have 'Set' in the name, since it sets the type of matrix to use.
1491fd1303e9SJungho Lee 
1492e69d4d44SBarry Smith     Level: intermediate
1493e69d4d44SBarry Smith 
1494fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC
1495e69d4d44SBarry Smith 
1496e69d4d44SBarry Smith @*/
14977087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1498e69d4d44SBarry Smith {
14994ac538c5SBarry Smith   PetscErrorCode ierr;
1500e69d4d44SBarry Smith 
1501e69d4d44SBarry Smith   PetscFunctionBegin;
15020700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
15034ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr);
1504e69d4d44SBarry Smith   PetscFunctionReturn(0);
1505e69d4d44SBarry Smith }
1506e69d4d44SBarry Smith 
1507e69d4d44SBarry Smith #undef __FUNCT__
1508470b340bSDmitry Karpeev #define __FUNCT__ "PCFieldSplitSchurGetS"
150945e7fc46SDmitry Karpeev /*@
1510470b340bSDmitry Karpeev     PCFieldSplitSchurGetS -  extract the MatSchurComplement object used by this PC in case it needs to be configured separately
151145e7fc46SDmitry Karpeev 
151245e7fc46SDmitry Karpeev     Not collective
151345e7fc46SDmitry Karpeev 
151445e7fc46SDmitry Karpeev     Input Parameter:
151545e7fc46SDmitry Karpeev .   pc      - the preconditioner context
151645e7fc46SDmitry Karpeev 
151745e7fc46SDmitry Karpeev     Output Parameter:
1518470b340bSDmitry Karpeev .   S       - the Schur complement matrix
151945e7fc46SDmitry Karpeev 
1520470b340bSDmitry Karpeev     Notes:
1521470b340bSDmitry Karpeev     This matrix should not be destroyed using MatDestroy(); rather, use PCFieldSplitSchurRestoreS().
152245e7fc46SDmitry Karpeev 
152345e7fc46SDmitry Karpeev     Level: advanced
152445e7fc46SDmitry Karpeev 
1525470b340bSDmitry Karpeev .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSchurPreType, PCFieldSplitSchurPrecondition(), MatSchurComplement, PCFieldSplitSchurRestoreS()
152645e7fc46SDmitry Karpeev 
152745e7fc46SDmitry Karpeev @*/
1528470b340bSDmitry Karpeev PetscErrorCode  PCFieldSplitSchurGetS(PC pc,Mat *S)
152945e7fc46SDmitry Karpeev {
153045e7fc46SDmitry Karpeev   PetscErrorCode ierr;
153145e7fc46SDmitry Karpeev   const char*    t;
153245e7fc46SDmitry Karpeev   PetscBool      isfs;
153345e7fc46SDmitry Karpeev   PC_FieldSplit  *jac;
153445e7fc46SDmitry Karpeev 
153545e7fc46SDmitry Karpeev   PetscFunctionBegin;
153645e7fc46SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
153745e7fc46SDmitry Karpeev   ierr = PetscObjectGetType((PetscObject)pc,&t);CHKERRQ(ierr);
153845e7fc46SDmitry Karpeev   ierr = PetscStrcmp(t,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
153945e7fc46SDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PC of type PCFIELDSPLIT, got %s instead",t);
154045e7fc46SDmitry Karpeev   jac = (PC_FieldSplit*)pc->data;
1541470b340bSDmitry 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);
1542470b340bSDmitry Karpeev   if (S) *S = jac->schur;
154345e7fc46SDmitry Karpeev   PetscFunctionReturn(0);
154445e7fc46SDmitry Karpeev }
154545e7fc46SDmitry Karpeev 
154645e7fc46SDmitry Karpeev #undef __FUNCT__
1547470b340bSDmitry Karpeev #define __FUNCT__ "PCFieldSplitSchurRestoreS"
1548470b340bSDmitry Karpeev /*@
1549470b340bSDmitry Karpeev     PCFieldSplitSchurRestoreS -  restores the MatSchurComplement object used by this PC
1550470b340bSDmitry Karpeev 
1551470b340bSDmitry Karpeev     Not collective
1552470b340bSDmitry Karpeev 
1553470b340bSDmitry Karpeev     Input Parameters:
1554470b340bSDmitry Karpeev +   pc      - the preconditioner context
1555470b340bSDmitry Karpeev .   S       - the Schur complement matrix
1556470b340bSDmitry Karpeev 
1557470b340bSDmitry Karpeev     Level: advanced
1558470b340bSDmitry Karpeev 
1559470b340bSDmitry Karpeev .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSchurPreType, PCFieldSplitSchurPrecondition(), MatSchurComplement, PCFieldSplitSchurGetS()
1560470b340bSDmitry Karpeev 
1561470b340bSDmitry Karpeev @*/
1562bca69d2bSDmitry Karpeev PetscErrorCode  PCFieldSplitSchurRestoreS(PC pc,Mat *S)
1563470b340bSDmitry Karpeev {
1564470b340bSDmitry Karpeev   PetscErrorCode ierr;
1565470b340bSDmitry Karpeev   const char*    t;
1566470b340bSDmitry Karpeev   PetscBool      isfs;
1567470b340bSDmitry Karpeev   PC_FieldSplit  *jac;
1568470b340bSDmitry Karpeev 
1569470b340bSDmitry Karpeev   PetscFunctionBegin;
1570470b340bSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1571470b340bSDmitry Karpeev   ierr = PetscObjectGetType((PetscObject)pc,&t);CHKERRQ(ierr);
1572470b340bSDmitry Karpeev   ierr = PetscStrcmp(t,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
1573470b340bSDmitry Karpeev   if (!isfs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected PC of type PCFIELDSPLIT, got %s instead",t);
1574470b340bSDmitry Karpeev   jac = (PC_FieldSplit*)pc->data;
1575470b340bSDmitry 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);
1576bca69d2bSDmitry Karpeev   if (!S || *S != jac->schur) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"MatSchurComplement restored is not the same as gotten");
1577470b340bSDmitry Karpeev   PetscFunctionReturn(0);
1578470b340bSDmitry Karpeev }
1579470b340bSDmitry Karpeev 
1580470b340bSDmitry Karpeev 
1581470b340bSDmitry Karpeev #undef __FUNCT__
1582e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
15831e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1584e69d4d44SBarry Smith {
1585e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1586084e4875SJed Brown   PetscErrorCode ierr;
1587e69d4d44SBarry Smith 
1588e69d4d44SBarry Smith   PetscFunctionBegin;
1589084e4875SJed Brown   jac->schurpre = ptype;
1590a7476a74SDmitry Karpeev   if (ptype == PC_FIELDSPLIT_SCHUR_PRE_USER && pre) {
15916bf464f9SBarry Smith     ierr            = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
1592084e4875SJed Brown     jac->schur_user = pre;
1593084e4875SJed Brown     ierr            = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1594084e4875SJed Brown   }
1595e69d4d44SBarry Smith   PetscFunctionReturn(0);
1596e69d4d44SBarry Smith }
1597e69d4d44SBarry Smith 
159830ad9308SMatthew Knepley #undef __FUNCT__
1599c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType"
1600ab1df9f5SJed Brown /*@
1601c9c6ffaaSJed Brown     PCFieldSplitSetSchurFactType -  sets which blocks of the approximate block factorization to retain
1602ab1df9f5SJed Brown 
1603ab1df9f5SJed Brown     Collective on PC
1604ab1df9f5SJed Brown 
1605ab1df9f5SJed Brown     Input Parameters:
1606ab1df9f5SJed Brown +   pc  - the preconditioner context
1607c9c6ffaaSJed Brown -   ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default
1608ab1df9f5SJed Brown 
1609ab1df9f5SJed Brown     Options Database:
1610c9c6ffaaSJed Brown .     -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full
1611ab1df9f5SJed Brown 
1612ab1df9f5SJed Brown 
1613ab1df9f5SJed Brown     Level: intermediate
1614ab1df9f5SJed Brown 
1615ab1df9f5SJed Brown     Notes:
1616ab1df9f5SJed Brown     The FULL factorization is
1617ab1df9f5SJed Brown 
1618ab1df9f5SJed Brown $   (A   B)  = (1       0) (A   0) (1  Ainv*B)
1619ab1df9f5SJed Brown $   (C   D)    (C*Ainv  1) (0   S) (0     1  )
1620ab1df9f5SJed Brown 
16216be4592eSBarry 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,
16226be4592eSBarry 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).
1623ab1df9f5SJed Brown 
16246be4592eSBarry 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
16256be4592eSBarry 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
16266be4592eSBarry 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,
16276be4592eSBarry 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.
1628ab1df9f5SJed Brown 
16296be4592eSBarry 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
16306be4592eSBarry 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).
1631ab1df9f5SJed Brown 
1632ab1df9f5SJed Brown     References:
1633ab1df9f5SJed Brown     Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000) pp. 1969-1972.
1634ab1df9f5SJed Brown 
1635ab1df9f5SJed Brown     Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001), pp. 1050-1051.
1636ab1df9f5SJed Brown 
1637ab1df9f5SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType
1638ab1df9f5SJed Brown @*/
1639c9c6ffaaSJed Brown PetscErrorCode  PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype)
1640ab1df9f5SJed Brown {
1641ab1df9f5SJed Brown   PetscErrorCode ierr;
1642ab1df9f5SJed Brown 
1643ab1df9f5SJed Brown   PetscFunctionBegin;
1644ab1df9f5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1645c9c6ffaaSJed Brown   ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr);
1646ab1df9f5SJed Brown   PetscFunctionReturn(0);
1647ab1df9f5SJed Brown }
1648ab1df9f5SJed Brown 
1649ab1df9f5SJed Brown #undef __FUNCT__
1650c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType_FieldSplit"
16511e6b0712SBarry Smith static PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype)
1652ab1df9f5SJed Brown {
1653ab1df9f5SJed Brown   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1654ab1df9f5SJed Brown 
1655ab1df9f5SJed Brown   PetscFunctionBegin;
1656ab1df9f5SJed Brown   jac->schurfactorization = ftype;
1657ab1df9f5SJed Brown   PetscFunctionReturn(0);
1658ab1df9f5SJed Brown }
1659ab1df9f5SJed Brown 
1660ab1df9f5SJed Brown #undef __FUNCT__
166130ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
166230ad9308SMatthew Knepley /*@C
16638c03b21aSDmitry Karpeev    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
166430ad9308SMatthew Knepley 
166530ad9308SMatthew Knepley    Collective on KSP
166630ad9308SMatthew Knepley 
166730ad9308SMatthew Knepley    Input Parameter:
166830ad9308SMatthew Knepley .  pc - the preconditioner context
166930ad9308SMatthew Knepley 
167030ad9308SMatthew Knepley    Output Parameters:
1671a04f6461SBarry Smith +  A00 - the (0,0) block
1672a04f6461SBarry Smith .  A01 - the (0,1) block
1673a04f6461SBarry Smith .  A10 - the (1,0) block
1674a04f6461SBarry Smith -  A11 - the (1,1) block
167530ad9308SMatthew Knepley 
167630ad9308SMatthew Knepley    Level: advanced
167730ad9308SMatthew Knepley 
167830ad9308SMatthew Knepley .seealso: PCFIELDSPLIT
167930ad9308SMatthew Knepley @*/
1680a04f6461SBarry Smith PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
168130ad9308SMatthew Knepley {
168230ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;
168330ad9308SMatthew Knepley 
168430ad9308SMatthew Knepley   PetscFunctionBegin;
16850700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1686ce94432eSBarry Smith   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
1687a04f6461SBarry Smith   if (A00) *A00 = jac->pmat[0];
1688a04f6461SBarry Smith   if (A01) *A01 = jac->B;
1689a04f6461SBarry Smith   if (A10) *A10 = jac->C;
1690a04f6461SBarry Smith   if (A11) *A11 = jac->pmat[1];
169130ad9308SMatthew Knepley   PetscFunctionReturn(0);
169230ad9308SMatthew Knepley }
169330ad9308SMatthew Knepley 
169479416396SBarry Smith #undef __FUNCT__
169579416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
16961e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
169779416396SBarry Smith {
169879416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1699e69d4d44SBarry Smith   PetscErrorCode ierr;
170079416396SBarry Smith 
170179416396SBarry Smith   PetscFunctionBegin;
170279416396SBarry Smith   jac->type = type;
17033b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
17043b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
17053b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
17062fa5cd67SKarl Rupp 
1707bdf89e91SBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1708bdf89e91SBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1709bdf89e91SBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr);
1710e69d4d44SBarry Smith 
17113b224e63SBarry Smith   } else {
17123b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
17133b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
17142fa5cd67SKarl Rupp 
1715bdf89e91SBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
1716bdf89e91SBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",0);CHKERRQ(ierr);
1717bdf89e91SBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetSchurFactType_C",0);CHKERRQ(ierr);
17183b224e63SBarry Smith   }
171979416396SBarry Smith   PetscFunctionReturn(0);
172079416396SBarry Smith }
172179416396SBarry Smith 
172251f519a2SBarry Smith #undef __FUNCT__
172351f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
17241e6b0712SBarry Smith static PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
172551f519a2SBarry Smith {
172651f519a2SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
172751f519a2SBarry Smith 
172851f519a2SBarry Smith   PetscFunctionBegin;
1729ce94432eSBarry Smith   if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
1730ce94432eSBarry 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);
173151f519a2SBarry Smith   jac->bs = bs;
173251f519a2SBarry Smith   PetscFunctionReturn(0);
173351f519a2SBarry Smith }
173451f519a2SBarry Smith 
173579416396SBarry Smith #undef __FUNCT__
173679416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
1737bc08b0f1SBarry Smith /*@
173879416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
173979416396SBarry Smith 
174079416396SBarry Smith    Collective on PC
174179416396SBarry Smith 
174279416396SBarry Smith    Input Parameter:
174379416396SBarry Smith .  pc - the preconditioner context
174481540f2fSBarry Smith .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
174579416396SBarry Smith 
174679416396SBarry Smith    Options Database Key:
1747a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
174879416396SBarry Smith 
1749b02e2d75SMatthew G Knepley    Level: Intermediate
175079416396SBarry Smith 
175179416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
175279416396SBarry Smith 
175379416396SBarry Smith .seealso: PCCompositeSetType()
175479416396SBarry Smith 
175579416396SBarry Smith @*/
17567087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
175779416396SBarry Smith {
17584ac538c5SBarry Smith   PetscErrorCode ierr;
175979416396SBarry Smith 
176079416396SBarry Smith   PetscFunctionBegin;
17610700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
17624ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr);
176379416396SBarry Smith   PetscFunctionReturn(0);
176479416396SBarry Smith }
176579416396SBarry Smith 
1766b02e2d75SMatthew G Knepley #undef __FUNCT__
1767b02e2d75SMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetType"
1768b02e2d75SMatthew G Knepley /*@
1769b02e2d75SMatthew G Knepley   PCFieldSplitGetType - Gets the type of fieldsplit preconditioner.
1770b02e2d75SMatthew G Knepley 
1771b02e2d75SMatthew G Knepley   Not collective
1772b02e2d75SMatthew G Knepley 
1773b02e2d75SMatthew G Knepley   Input Parameter:
1774b02e2d75SMatthew G Knepley . pc - the preconditioner context
1775b02e2d75SMatthew G Knepley 
1776b02e2d75SMatthew G Knepley   Output Parameter:
1777b02e2d75SMatthew G Knepley . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
1778b02e2d75SMatthew G Knepley 
1779b02e2d75SMatthew G Knepley   Level: Intermediate
1780b02e2d75SMatthew G Knepley 
1781b02e2d75SMatthew G Knepley .keywords: PC, set, type, composite preconditioner, additive, multiplicative
1782b02e2d75SMatthew G Knepley .seealso: PCCompositeSetType()
1783b02e2d75SMatthew G Knepley @*/
1784b02e2d75SMatthew G Knepley PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type)
1785b02e2d75SMatthew G Knepley {
1786b02e2d75SMatthew G Knepley   PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;
1787b02e2d75SMatthew G Knepley 
1788b02e2d75SMatthew G Knepley   PetscFunctionBegin;
1789b02e2d75SMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1790b02e2d75SMatthew G Knepley   PetscValidIntPointer(type,2);
1791b02e2d75SMatthew G Knepley   *type = jac->type;
1792b02e2d75SMatthew G Knepley   PetscFunctionReturn(0);
1793b02e2d75SMatthew G Knepley }
1794b02e2d75SMatthew G Knepley 
17954ab8060aSDmitry Karpeev #undef __FUNCT__
17964ab8060aSDmitry Karpeev #define __FUNCT__ "PCFieldSplitSetDMSplits"
17974ab8060aSDmitry Karpeev /*@
17984ab8060aSDmitry Karpeev    PCFieldSplitSetDMSplits - Flags whether DMCreateFieldDecomposition() should be used to define the splits, whenever possible.
17994ab8060aSDmitry Karpeev 
18004ab8060aSDmitry Karpeev    Logically Collective
18014ab8060aSDmitry Karpeev 
18024ab8060aSDmitry Karpeev    Input Parameters:
18034ab8060aSDmitry Karpeev +  pc   - the preconditioner context
18044ab8060aSDmitry Karpeev -  flg  - boolean indicating whether to use field splits defined by the DM
18054ab8060aSDmitry Karpeev 
18064ab8060aSDmitry Karpeev    Options Database Key:
18074ab8060aSDmitry Karpeev .  -pc_fieldsplit_dm_splits
18084ab8060aSDmitry Karpeev 
18094ab8060aSDmitry Karpeev    Level: Intermediate
18104ab8060aSDmitry Karpeev 
18114ab8060aSDmitry Karpeev .keywords: PC, DM, composite preconditioner, additive, multiplicative
18124ab8060aSDmitry Karpeev 
18134ab8060aSDmitry Karpeev .seealso: PCFieldSplitGetDMSplits()
18144ab8060aSDmitry Karpeev 
18154ab8060aSDmitry Karpeev @*/
18164ab8060aSDmitry Karpeev PetscErrorCode  PCFieldSplitSetDMSplits(PC pc,PetscBool flg)
18174ab8060aSDmitry Karpeev {
18184ab8060aSDmitry Karpeev   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
18194ab8060aSDmitry Karpeev   PetscBool      isfs;
18204ab8060aSDmitry Karpeev   PetscErrorCode ierr;
18214ab8060aSDmitry Karpeev 
18224ab8060aSDmitry Karpeev   PetscFunctionBegin;
18234ab8060aSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
18244ab8060aSDmitry Karpeev   PetscValidLogicalCollectiveBool(pc,flg,2);
18254ab8060aSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
18264ab8060aSDmitry Karpeev   if (isfs) {
18274ab8060aSDmitry Karpeev     jac->dm_splits = flg;
18284ab8060aSDmitry Karpeev   }
18294ab8060aSDmitry Karpeev   PetscFunctionReturn(0);
18304ab8060aSDmitry Karpeev }
18314ab8060aSDmitry Karpeev 
18324ab8060aSDmitry Karpeev 
18334ab8060aSDmitry Karpeev #undef __FUNCT__
18344ab8060aSDmitry Karpeev #define __FUNCT__ "PCFieldSplitGetDMSplits"
18354ab8060aSDmitry Karpeev /*@
18364ab8060aSDmitry Karpeev    PCFieldSplitGetDMSplits - Returns flag indicating whether DMCreateFieldDecomposition() should be used to define the splits, whenever possible.
18374ab8060aSDmitry Karpeev 
18384ab8060aSDmitry Karpeev    Logically Collective
18394ab8060aSDmitry Karpeev 
18404ab8060aSDmitry Karpeev    Input Parameter:
18414ab8060aSDmitry Karpeev .  pc   - the preconditioner context
18424ab8060aSDmitry Karpeev 
18434ab8060aSDmitry Karpeev    Output Parameter:
18444ab8060aSDmitry Karpeev .  flg  - boolean indicating whether to use field splits defined by the DM
18454ab8060aSDmitry Karpeev 
18464ab8060aSDmitry Karpeev    Level: Intermediate
18474ab8060aSDmitry Karpeev 
18484ab8060aSDmitry Karpeev .keywords: PC, DM, composite preconditioner, additive, multiplicative
18494ab8060aSDmitry Karpeev 
18504ab8060aSDmitry Karpeev .seealso: PCFieldSplitSetDMSplits()
18514ab8060aSDmitry Karpeev 
18524ab8060aSDmitry Karpeev @*/
18534ab8060aSDmitry Karpeev PetscErrorCode  PCFieldSplitGetDMSplits(PC pc,PetscBool* flg)
18544ab8060aSDmitry Karpeev {
18554ab8060aSDmitry Karpeev   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
18564ab8060aSDmitry Karpeev   PetscBool      isfs;
18574ab8060aSDmitry Karpeev   PetscErrorCode ierr;
18584ab8060aSDmitry Karpeev 
18594ab8060aSDmitry Karpeev   PetscFunctionBegin;
18604ab8060aSDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
18614ab8060aSDmitry Karpeev   PetscValidPointer(flg,2);
18624ab8060aSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCFIELDSPLIT,&isfs);CHKERRQ(ierr);
18634ab8060aSDmitry Karpeev   if (isfs) {
18644ab8060aSDmitry Karpeev     if(flg) *flg = jac->dm_splits;
18654ab8060aSDmitry Karpeev   }
18664ab8060aSDmitry Karpeev   PetscFunctionReturn(0);
18674ab8060aSDmitry Karpeev }
18684ab8060aSDmitry Karpeev 
18690971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
18700971522cSBarry Smith /*MC
1871a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
1872a04f6461SBarry Smith                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.
18730971522cSBarry Smith 
1874edf189efSBarry Smith      To set options on the solvers for each block append -fieldsplit_ to all the PC
1875edf189efSBarry Smith         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
18760971522cSBarry Smith 
1877a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
187869a612a9SBarry Smith          and set the options directly on the resulting KSP object
18790971522cSBarry Smith 
18800971522cSBarry Smith    Level: intermediate
18810971522cSBarry Smith 
188279416396SBarry Smith    Options Database Keys:
188381540f2fSBarry Smith +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
188481540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
188581540f2fSBarry Smith                               been supplied explicitly by -pc_fieldsplit_%d_fields
188681540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
18870f188ba9SJed Brown .   -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur> - type of relaxation or factorization splitting
1888*2e71c61dSMatthew G. Knepley .   -pc_fieldsplit_schur_precondition <self,selfp,user,a11,full> - default is a11
1889435f959eSBarry Smith .   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver
189079416396SBarry Smith 
18915d4c12cdSJungho Lee -    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
18925d4c12cdSJungho Lee      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
18935d4c12cdSJungho Lee 
1894c8a0d604SMatthew G Knepley    Notes:
1895c8a0d604SMatthew G Knepley     Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1896d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
1897d32f9abdSBarry Smith 
1898d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
1899d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1900d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
1901d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
1902d32f9abdSBarry Smith 
1903c8a0d604SMatthew G Knepley $     For the Schur complement preconditioner if J = ( A00 A01 )
1904c8a0d604SMatthew G Knepley $                                                    ( A10 A11 )
1905c8a0d604SMatthew G Knepley $     the preconditioner using full factorization is
190613b05affSJed Brown $              ( I   -ksp(A00) A01 ) ( inv(A00)     0  ) (     I          0  )
1907c8a0d604SMatthew G Knepley $              ( 0         I       ) (   0      ksp(S) ) ( -A10 ksp(A00)  I  )
190813b05affSJed Brown      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_.  S is the Schur complement
190913b05affSJed Brown $              S = A11 - A10 ksp(A00) A01
191013b05affSJed 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
1911c8a0d604SMatthew G Knepley      in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0,
1912c8a0d604SMatthew G Knepley      it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
1913a04f6461SBarry Smith      A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
1914a7476a74SDmitry Karpeev      option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc.
1915a7476a74SDmitry Karpeev      When option -fieldsplit_schur_precondition selfp is given, an approximation to S is assembled --
1916a7476a74SDmitry Karpeev      Sp = A11 - A10 inv(diag(A00)) A01, which has type AIJ and can be used with a variety of preconditioners
1917a7476a74SDmitry Karpeev      (e.g., -fieldsplit_1_pc_type asm).
1918a7476a74SDmitry Karpeev      The factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above,
19195668aaf4SBarry Smith      diag gives
1920c8a0d604SMatthew G Knepley $              ( inv(A00)     0   )
1921c8a0d604SMatthew G Knepley $              (   0      -ksp(S) )
19225668aaf4SBarry 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
1923c8a0d604SMatthew G Knepley $              (  A00   0 )
1924c8a0d604SMatthew G Knepley $              (  A10   S )
1925c8a0d604SMatthew G Knepley      where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of
1926c8a0d604SMatthew G Knepley $              ( A00 A01 )
1927c8a0d604SMatthew G Knepley $              (  0   S  )
1928c8a0d604SMatthew G Knepley      where again the inverses of A00 and S are applied using KSPs.
1929e69d4d44SBarry Smith 
1930edf189efSBarry Smith      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
1931edf189efSBarry Smith      is used automatically for a second block.
1932edf189efSBarry Smith 
1933ff218e97SBarry Smith      The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1.
1934ff218e97SBarry Smith      Generally it should be used with the AIJ format.
1935ff218e97SBarry Smith 
1936ff218e97SBarry Smith      The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see,
1937ff218e97SBarry Smith      for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT
1938ff218e97SBarry Smith      inside a smoother resulting in "Distributive Smoothers".
19390716a85fSBarry Smith 
1940a541d17aSBarry Smith    Concepts: physics based preconditioners, block preconditioners
19410971522cSBarry Smith 
19427e8cb189SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
1943e69d4d44SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
19440971522cSBarry Smith M*/
19450971522cSBarry Smith 
19460971522cSBarry Smith #undef __FUNCT__
19470971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
19488cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_FieldSplit(PC pc)
19490971522cSBarry Smith {
19500971522cSBarry Smith   PetscErrorCode ierr;
19510971522cSBarry Smith   PC_FieldSplit  *jac;
19520971522cSBarry Smith 
19530971522cSBarry Smith   PetscFunctionBegin;
1954b00a9115SJed Brown   ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr);
19552fa5cd67SKarl Rupp 
19560971522cSBarry Smith   jac->bs                 = -1;
19570971522cSBarry Smith   jac->nsplits            = 0;
19583e197d65SBarry Smith   jac->type               = PC_COMPOSITE_MULTIPLICATIVE;
1959e6cab6aaSJed Brown   jac->schurpre           = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
1960c9c6ffaaSJed Brown   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL;
1961fbe7908bSJed Brown   jac->dm_splits          = PETSC_TRUE;
196251f519a2SBarry Smith 
19630971522cSBarry Smith   pc->data = (void*)jac;
19640971522cSBarry Smith 
19650971522cSBarry Smith   pc->ops->apply           = PCApply_FieldSplit;
1966421e10b8SBarry Smith   pc->ops->applytranspose  = PCApplyTranspose_FieldSplit;
19670971522cSBarry Smith   pc->ops->setup           = PCSetUp_FieldSplit;
1968574deadeSBarry Smith   pc->ops->reset           = PCReset_FieldSplit;
19690971522cSBarry Smith   pc->ops->destroy         = PCDestroy_FieldSplit;
19700971522cSBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_FieldSplit;
19710971522cSBarry Smith   pc->ops->view            = PCView_FieldSplit;
19720971522cSBarry Smith   pc->ops->applyrichardson = 0;
19730971522cSBarry Smith 
1974bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
1975bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetFields_C",PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1976bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetIS_C",PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
1977bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetType_C",PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
1978bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
19790971522cSBarry Smith   PetscFunctionReturn(0);
19800971522cSBarry Smith }
19810971522cSBarry Smith 
19820971522cSBarry Smith 
1983a541d17aSBarry Smith 
1984