xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision ce94432eddcd14845bc7e8083b7f8ea723b9bf7d)
1dba47a55SKris Buschelman 
2b45d2f2cSJed Brown #include <petsc-private/pcimpl.h>     /*I "petscpc.h" I*/
307475bc1SBarry Smith #include <petscdmcomposite.h>
40971522cSBarry Smith 
5443836d0SMatthew G Knepley /*
6443836d0SMatthew G Knepley   There is a nice discussion of block preconditioners in
7443836d0SMatthew G Knepley 
8443836d0SMatthew G Knepley [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 
13e87fab1bSBarry Smith const char *const PCFieldSplitSchurPreTypes[] = {"SELF","A11","USER","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 */
32ace3abfcSBarry Smith   PetscBool       realdiagonal;                    /* Flag to use the diagonal blocks of mat preconditioned by pmat, instead of just pmat */
3330ad9308SMatthew Knepley   PetscInt        bs;                              /* Block size for IS and Mat structures */
3430ad9308SMatthew Knepley   PetscInt        nsplits;                         /* Number of field divisions defined */
3579416396SBarry Smith   Vec             *x,*y,w1,w2;
36519d70e2SJed Brown   Mat             *mat;                            /* The diagonal block for each split */
37519d70e2SJed Brown   Mat             *pmat;                           /* The preconditioning diagonal block for each split */
3830ad9308SMatthew Knepley   Mat             *Afield;                         /* The rows of the matrix associated with each split */
39ace3abfcSBarry Smith   PetscBool       issetup;
402fa5cd67SKarl Rupp 
4130ad9308SMatthew Knepley   /* Only used when Schur complement preconditioning is used */
4230ad9308SMatthew Knepley   Mat                       B;                     /* The (0,1) block */
4330ad9308SMatthew Knepley   Mat                       C;                     /* The (1,0) block */
44443836d0SMatthew G Knepley   Mat                       schur;                 /* The Schur complement S = A11 - A10 A00^{-1} A01, the KSP here, kspinner, is H_1 in [El08] */
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 */
530971522cSBarry Smith } PC_FieldSplit;
540971522cSBarry Smith 
5516913363SBarry Smith /*
5616913363SBarry Smith     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
5716913363SBarry Smith    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
5816913363SBarry Smith    PC you could change this.
5916913363SBarry Smith */
60084e4875SJed Brown 
61e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it.  This way the
62084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */
63084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac)
64084e4875SJed Brown {
65084e4875SJed Brown   switch (jac->schurpre) {
66084e4875SJed Brown   case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur;
67e87fab1bSBarry Smith   case PC_FIELDSPLIT_SCHUR_PRE_A11: return jac->pmat[1];
68084e4875SJed Brown   case PC_FIELDSPLIT_SCHUR_PRE_USER:   /* Use a user-provided matrix if it is given, otherwise diagonal block */
69084e4875SJed Brown   default:
70084e4875SJed Brown     return jac->schur_user ? jac->schur_user : jac->pmat[1];
71084e4875SJed Brown   }
72084e4875SJed Brown }
73084e4875SJed Brown 
74084e4875SJed Brown 
750971522cSBarry Smith #undef __FUNCT__
760971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
770971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
780971522cSBarry Smith {
790971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
800971522cSBarry Smith   PetscErrorCode    ierr;
81d9884438SBarry Smith   PetscBool         iascii,isdraw;
820971522cSBarry Smith   PetscInt          i,j;
835a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
840971522cSBarry Smith 
850971522cSBarry Smith   PetscFunctionBegin;
86251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
87d9884438SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
880971522cSBarry Smith   if (iascii) {
892c9966d7SBarry Smith     if (jac->bs > 0) {
9051f519a2SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
912c9966d7SBarry Smith     } else {
922c9966d7SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D\n",PCCompositeTypes[jac->type],jac->nsplits);CHKERRQ(ierr);
932c9966d7SBarry Smith     }
94a3df900dSMatthew G Knepley     if (jac->realdiagonal) {
95a3df900dSMatthew G Knepley       ierr = PetscViewerASCIIPrintf(viewer,"  using actual matrix for blocks rather than preconditioner matrix\n");CHKERRQ(ierr);
96a3df900dSMatthew G Knepley     }
9769a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
980971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
990971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
1001ab39975SBarry Smith       if (ilink->fields) {
1010971522cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
10279416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
1035a9f2f41SSatish Balay         for (j=0; j<ilink->nfields; j++) {
10479416396SBarry Smith           if (j > 0) {
10579416396SBarry Smith             ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
10679416396SBarry Smith           }
1075a9f2f41SSatish Balay           ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
1080971522cSBarry Smith         }
1090971522cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
11079416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1111ab39975SBarry Smith       } else {
1121ab39975SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1131ab39975SBarry Smith       }
1145a9f2f41SSatish Balay       ierr  = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
1155a9f2f41SSatish Balay       ilink = ilink->next;
1160971522cSBarry Smith     }
1170971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1182fa5cd67SKarl Rupp   }
1192fa5cd67SKarl Rupp 
1202fa5cd67SKarl Rupp  if (isdraw) {
121d9884438SBarry Smith     PetscDraw draw;
122d9884438SBarry Smith     PetscReal x,y,w,wd;
123d9884438SBarry Smith 
124d9884438SBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
125d9884438SBarry Smith     ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
126d9884438SBarry Smith     w    = 2*PetscMin(1.0 - x,x);
127d9884438SBarry Smith     wd   = w/(jac->nsplits + 1);
128d9884438SBarry Smith     x    = x - wd*(jac->nsplits-1)/2.0;
129d9884438SBarry Smith     for (i=0; i<jac->nsplits; i++) {
130d9884438SBarry Smith       ierr  = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
131d9884438SBarry Smith       ierr  = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
132d9884438SBarry Smith       ierr  = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
133d9884438SBarry Smith       x    += wd;
134d9884438SBarry Smith       ilink = ilink->next;
135d9884438SBarry Smith     }
1360971522cSBarry Smith   }
1370971522cSBarry Smith   PetscFunctionReturn(0);
1380971522cSBarry Smith }
1390971522cSBarry Smith 
1400971522cSBarry Smith #undef __FUNCT__
1413b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur"
1423b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
1433b224e63SBarry Smith {
1443b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1453b224e63SBarry Smith   PetscErrorCode    ierr;
1464996c5bdSBarry Smith   PetscBool         iascii,isdraw;
1473b224e63SBarry Smith   PetscInt          i,j;
1483b224e63SBarry Smith   PC_FieldSplitLink ilink = jac->head;
1493b224e63SBarry Smith 
1503b224e63SBarry Smith   PetscFunctionBegin;
151251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1524996c5bdSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
1533b224e63SBarry Smith   if (iascii) {
1542c9966d7SBarry Smith     if (jac->bs > 0) {
155c9c6ffaaSJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
1562c9966d7SBarry Smith     } else {
1572c9966d7SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, factorization %s\n",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
1582c9966d7SBarry Smith     }
159a3df900dSMatthew G Knepley     if (jac->realdiagonal) {
160a3df900dSMatthew G Knepley       ierr = PetscViewerASCIIPrintf(viewer,"  using actual matrix for blocks rather than preconditioner matrix\n");CHKERRQ(ierr);
161a3df900dSMatthew G Knepley     }
1623e8b8b31SMatthew G Knepley     switch (jac->schurpre) {
1633e8b8b31SMatthew G Knepley     case PC_FIELDSPLIT_SCHUR_PRE_SELF:
1643e8b8b31SMatthew G Knepley       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from S itself\n");CHKERRQ(ierr);break;
165e87fab1bSBarry Smith     case PC_FIELDSPLIT_SCHUR_PRE_A11:
166e87fab1bSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from A11\n");CHKERRQ(ierr);break;
1673e8b8b31SMatthew G Knepley     case PC_FIELDSPLIT_SCHUR_PRE_USER:
1683e8b8b31SMatthew G Knepley       if (jac->schur_user) {
1693e8b8b31SMatthew G Knepley         ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from user provided matrix\n");CHKERRQ(ierr);
1703e8b8b31SMatthew G Knepley       } else {
171e87fab1bSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioner for the Schur complement formed from A11\n");CHKERRQ(ierr);
1723e8b8b31SMatthew G Knepley       }
1733e8b8b31SMatthew G Knepley       break;
1743e8b8b31SMatthew G Knepley     default:
1753e8b8b31SMatthew G Knepley       SETERRQ1(((PetscObject) pc)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid Schur preconditioning type: %d", jac->schurpre);
1763e8b8b31SMatthew G Knepley     }
1773b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Split info:\n");CHKERRQ(ierr);
1783b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1793b224e63SBarry Smith     for (i=0; i<jac->nsplits; i++) {
1803b224e63SBarry Smith       if (ilink->fields) {
1813b224e63SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
1823b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
1833b224e63SBarry Smith         for (j=0; j<ilink->nfields; j++) {
1843b224e63SBarry Smith           if (j > 0) {
1853b224e63SBarry Smith             ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
1863b224e63SBarry Smith           }
1873b224e63SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
1883b224e63SBarry Smith         }
1893b224e63SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1903b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1913b224e63SBarry Smith       } else {
1923b224e63SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1933b224e63SBarry Smith       }
1943b224e63SBarry Smith       ilink = ilink->next;
1953b224e63SBarry Smith     }
196435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block \n");CHKERRQ(ierr);
1973b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
198443836d0SMatthew G Knepley     ierr = KSPView(jac->head->ksp,viewer);CHKERRQ(ierr);
1993b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
200443836d0SMatthew G Knepley     if (jac->kspupper != jac->head->ksp) {
201443836d0SMatthew G Knepley       ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for upper A00 in upper triangular factor \n");CHKERRQ(ierr);
202443836d0SMatthew G Knepley       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
203443836d0SMatthew G Knepley       ierr = KSPView(jac->kspupper,viewer);CHKERRQ(ierr);
204443836d0SMatthew G Knepley       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
205443836d0SMatthew G Knepley     }
206435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr);
2073b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
20812cae6f2SJed Brown     if (jac->kspschur) {
2093b224e63SBarry Smith       ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
21012cae6f2SJed Brown     } else {
21112cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
21212cae6f2SJed Brown     }
2133b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
2143b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
2154996c5bdSBarry Smith   } else if (isdraw) {
2164996c5bdSBarry Smith     PetscDraw draw;
2174996c5bdSBarry Smith     PetscReal x,y,w,wd,h;
2184996c5bdSBarry Smith     PetscInt  cnt = 2;
2194996c5bdSBarry Smith     char      str[32];
2204996c5bdSBarry Smith 
2214996c5bdSBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
2224996c5bdSBarry Smith     ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
223c74581afSBarry Smith     if (jac->kspupper != jac->head->ksp) cnt++;
224c74581afSBarry Smith     w  = 2*PetscMin(1.0 - x,x);
225c74581afSBarry Smith     wd = w/(cnt + 1);
226c74581afSBarry Smith 
2274996c5bdSBarry Smith     ierr = PetscSNPrintf(str,32,"Schur fact. %s",PCFieldSplitSchurFactTypes[jac->schurfactorization]);CHKERRQ(ierr);
2280298fd71SBarry Smith     ierr = PetscDrawBoxedString(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr);
2294996c5bdSBarry Smith     y   -= h;
2304996c5bdSBarry Smith     if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_USER &&  !jac->schur_user) {
231e87fab1bSBarry Smith       ierr = PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[PC_FIELDSPLIT_SCHUR_PRE_A11]);CHKERRQ(ierr);
2323b224e63SBarry Smith     } else {
2334996c5bdSBarry Smith       ierr = PetscSNPrintf(str,32,"Prec. for Schur from %s",PCFieldSplitSchurPreTypes[jac->schurpre]);CHKERRQ(ierr);
2344996c5bdSBarry Smith     }
2350298fd71SBarry Smith     ierr = PetscDrawBoxedString(draw,x+wd*(cnt-1)/2.0,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr);
2364996c5bdSBarry Smith     y   -= h;
2374996c5bdSBarry Smith     x    = x - wd*(cnt-1)/2.0;
2384996c5bdSBarry Smith 
2394996c5bdSBarry Smith     ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
2404996c5bdSBarry Smith     ierr = KSPView(jac->head->ksp,viewer);CHKERRQ(ierr);
2414996c5bdSBarry Smith     ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
2424996c5bdSBarry Smith     if (jac->kspupper != jac->head->ksp) {
2434996c5bdSBarry Smith       x   += wd;
2444996c5bdSBarry Smith       ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
2454996c5bdSBarry Smith       ierr = KSPView(jac->kspupper,viewer);CHKERRQ(ierr);
2464996c5bdSBarry Smith       ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
2474996c5bdSBarry Smith     }
2484996c5bdSBarry Smith     x   += wd;
2494996c5bdSBarry Smith     ierr = PetscDrawPushCurrentPoint(draw,x,y);CHKERRQ(ierr);
2504996c5bdSBarry Smith     ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
2514996c5bdSBarry Smith     ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
2523b224e63SBarry Smith   }
2533b224e63SBarry Smith   PetscFunctionReturn(0);
2543b224e63SBarry Smith }
2553b224e63SBarry Smith 
2563b224e63SBarry Smith #undef __FUNCT__
2576c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private"
2586c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */
2596c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc)
2606c924f48SJed Brown {
2616c924f48SJed Brown   PetscErrorCode ierr;
2626c924f48SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
2635d4c12cdSJungho Lee   PetscInt       i,nfields,*ifields,nfields_col,*ifields_col;
2645d4c12cdSJungho Lee   PetscBool      flg,flg_col;
2655d4c12cdSJungho Lee   char           optionname[128],splitname[8],optionname_col[128];
2666c924f48SJed Brown 
2676c924f48SJed Brown   PetscFunctionBegin;
2686c924f48SJed Brown   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr);
2695d4c12cdSJungho Lee   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields_col);CHKERRQ(ierr);
2706c924f48SJed Brown   for (i=0,flg=PETSC_TRUE;; i++) {
2718caf3d72SBarry Smith     ierr        = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr);
2728caf3d72SBarry Smith     ierr        = PetscSNPrintf(optionname,sizeof(optionname),"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr);
2738caf3d72SBarry Smith     ierr        = PetscSNPrintf(optionname_col,sizeof(optionname_col),"-pc_fieldsplit_%D_fields_col",i);CHKERRQ(ierr);
2746c924f48SJed Brown     nfields     = jac->bs;
27529499fbbSJungho Lee     nfields_col = jac->bs;
2766c924f48SJed Brown     ierr        = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr);
2775d4c12cdSJungho Lee     ierr        = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);CHKERRQ(ierr);
2786c924f48SJed Brown     if (!flg) break;
2795d4c12cdSJungho Lee     else if (flg && !flg_col) {
2805d4c12cdSJungho Lee       if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
2815d4c12cdSJungho Lee       ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);CHKERRQ(ierr);
2822fa5cd67SKarl Rupp     } else {
2835d4c12cdSJungho Lee       if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
2845d4c12cdSJungho Lee       if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match");
2855d4c12cdSJungho Lee       ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);CHKERRQ(ierr);
2865d4c12cdSJungho Lee     }
2876c924f48SJed Brown   }
2886c924f48SJed Brown   if (i > 0) {
2896c924f48SJed Brown     /* Makes command-line setting of splits take precedence over setting them in code.
2906c924f48SJed Brown        Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would
2916c924f48SJed Brown        create new splits, which would probably not be what the user wanted. */
2926c924f48SJed Brown     jac->splitdefined = PETSC_TRUE;
2936c924f48SJed Brown   }
2946c924f48SJed Brown   ierr = PetscFree(ifields);CHKERRQ(ierr);
2955d4c12cdSJungho Lee   ierr = PetscFree(ifields_col);CHKERRQ(ierr);
2966c924f48SJed Brown   PetscFunctionReturn(0);
2976c924f48SJed Brown }
2986c924f48SJed Brown 
2996c924f48SJed Brown #undef __FUNCT__
30069a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
30169a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
3020971522cSBarry Smith {
3030971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
3040971522cSBarry Smith   PetscErrorCode    ierr;
3055a9f2f41SSatish Balay   PC_FieldSplitLink ilink              = jac->head;
3067287d2a3SDmitry Karpeev   PetscBool         fieldsplit_default = PETSC_FALSE,stokes = PETSC_FALSE;
3076c924f48SJed Brown   PetscInt          i;
3080971522cSBarry Smith 
3090971522cSBarry Smith   PetscFunctionBegin;
3107287d2a3SDmitry Karpeev   /*
3117287d2a3SDmitry Karpeev    Kinda messy, but at least this now uses DMCreateFieldDecomposition() even with jac->reset.
3127287d2a3SDmitry Karpeev    Should probably be rewritten.
3137287d2a3SDmitry Karpeev    */
3147287d2a3SDmitry Karpeev   if (!ilink || jac->reset) {
3150298fd71SBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);CHKERRQ(ierr);
316d0af7cd3SBarry Smith     if (pc->dm && !stokes) {
317bafc1b83SMatthew G Knepley       PetscInt  numFields, f, i, j;
3180784a22cSJed Brown       char      **fieldNames;
3197b62db95SJungho Lee       IS        *fields;
320e7c4fc90SDmitry Karpeev       DM        *dms;
321bafc1b83SMatthew G Knepley       DM        subdm[128];
322bafc1b83SMatthew G Knepley       PetscBool flg;
323bafc1b83SMatthew G Knepley 
32416621825SDmitry Karpeev       ierr = DMCreateFieldDecomposition(pc->dm, &numFields, &fieldNames, &fields, &dms);CHKERRQ(ierr);
325bafc1b83SMatthew G Knepley       /* Allow the user to prescribe the splits */
326bafc1b83SMatthew G Knepley       for (i = 0, flg = PETSC_TRUE;; i++) {
327bafc1b83SMatthew G Knepley         PetscInt ifields[128];
328bafc1b83SMatthew G Knepley         IS       compField;
329bafc1b83SMatthew G Knepley         char     optionname[128], splitname[8];
330bafc1b83SMatthew G Knepley         PetscInt nfields = numFields;
331bafc1b83SMatthew G Knepley 
3328caf3d72SBarry Smith         ierr = PetscSNPrintf(optionname, sizeof(optionname), "-pc_fieldsplit_%D_fields", i);CHKERRQ(ierr);
333bafc1b83SMatthew G Knepley         ierr = PetscOptionsGetIntArray(((PetscObject) pc)->prefix, optionname, ifields, &nfields, &flg);CHKERRQ(ierr);
334bafc1b83SMatthew G Knepley         if (!flg) break;
335bafc1b83SMatthew G Knepley         if (numFields > 128) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Cannot currently support %d > 128 fields", numFields);
336bafc1b83SMatthew G Knepley         ierr = DMCreateSubDM(pc->dm, nfields, ifields, &compField, &subdm[i]);CHKERRQ(ierr);
337bafc1b83SMatthew G Knepley         if (nfields == 1) {
338bafc1b83SMatthew G Knepley           ierr = PCFieldSplitSetIS(pc, fieldNames[ifields[0]], compField);CHKERRQ(ierr);
3394d2389d7SMatthew G Knepley           /* ierr = PetscPrintf(((PetscObject) pc)->comm, "%s Field Indices:", fieldNames[ifields[0]]);CHKERRQ(ierr);
3400298fd71SBarry Smith              ierr = ISView(compField, NULL);CHKERRQ(ierr); */
341bafc1b83SMatthew G Knepley         } else {
3428caf3d72SBarry Smith           ierr = PetscSNPrintf(splitname, sizeof(splitname), "%D", i);CHKERRQ(ierr);
343bafc1b83SMatthew G Knepley           ierr = PCFieldSplitSetIS(pc, splitname, compField);CHKERRQ(ierr);
3444d2389d7SMatthew G Knepley           /* ierr = PetscPrintf(((PetscObject) pc)->comm, "%s Field Indices:", splitname);CHKERRQ(ierr);
3450298fd71SBarry Smith              ierr = ISView(compField, NULL);CHKERRQ(ierr); */
3467287d2a3SDmitry Karpeev         }
347bafc1b83SMatthew G Knepley         ierr = ISDestroy(&compField);CHKERRQ(ierr);
348bafc1b83SMatthew G Knepley         for (j = 0; j < nfields; ++j) {
349bafc1b83SMatthew G Knepley           f    = ifields[j];
3507b62db95SJungho Lee           ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr);
3517b62db95SJungho Lee           ierr = ISDestroy(&fields[f]);CHKERRQ(ierr);
3527b62db95SJungho Lee         }
353bafc1b83SMatthew G Knepley       }
354bafc1b83SMatthew G Knepley       if (i == 0) {
355bafc1b83SMatthew G Knepley         for (f = 0; f < numFields; ++f) {
356bafc1b83SMatthew G Knepley           ierr = PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);CHKERRQ(ierr);
357bafc1b83SMatthew G Knepley           ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr);
358bafc1b83SMatthew G Knepley           ierr = ISDestroy(&fields[f]);CHKERRQ(ierr);
359bafc1b83SMatthew G Knepley         }
360bafc1b83SMatthew G Knepley       } else {
361d724dfffSBarry Smith         for (j=0; j<numFields; j++) {
362d724dfffSBarry Smith           ierr = DMDestroy(dms+j);CHKERRQ(ierr);
363d724dfffSBarry Smith         }
364d724dfffSBarry Smith         ierr = PetscFree(dms);CHKERRQ(ierr);
365bafc1b83SMatthew G Knepley         ierr = PetscMalloc(i * sizeof(DM), &dms);CHKERRQ(ierr);
3662fa5cd67SKarl Rupp         for (j = 0; j < i; ++j) dms[j] = subdm[j];
367bafc1b83SMatthew G Knepley       }
3687b62db95SJungho Lee       ierr = PetscFree(fieldNames);CHKERRQ(ierr);
3697b62db95SJungho Lee       ierr = PetscFree(fields);CHKERRQ(ierr);
370e7c4fc90SDmitry Karpeev       if (dms) {
3718b8307b2SJed Brown         ierr = PetscInfo(pc, "Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
372bafc1b83SMatthew G Knepley         for (ilink = jac->head, i = 0; ilink; ilink = ilink->next, ++i) {
3737287d2a3SDmitry Karpeev           const char *prefix;
3747287d2a3SDmitry Karpeev           ierr = PetscObjectGetOptionsPrefix((PetscObject)(ilink->ksp),&prefix);CHKERRQ(ierr);
3757287d2a3SDmitry Karpeev           ierr = PetscObjectSetOptionsPrefix((PetscObject)(dms[i]), prefix);CHKERRQ(ierr);
3767b62db95SJungho Lee           ierr = KSPSetDM(ilink->ksp, dms[i]);CHKERRQ(ierr);
3777b62db95SJungho Lee           ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE);CHKERRQ(ierr);
3787287d2a3SDmitry Karpeev           ierr = PetscObjectIncrementTabLevel((PetscObject)dms[i],(PetscObject)ilink->ksp,0);CHKERRQ(ierr);
379e7c4fc90SDmitry Karpeev           ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);
3802fa5ba8aSJed Brown         }
3817b62db95SJungho Lee         ierr = PetscFree(dms);CHKERRQ(ierr);
3828b8307b2SJed Brown       }
38366ffff09SJed Brown     } else {
384521d7252SBarry Smith       if (jac->bs <= 0) {
385704ba839SBarry Smith         if (pc->pmat) {
386521d7252SBarry Smith           ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
3872fa5cd67SKarl Rupp         } else jac->bs = 1;
388521d7252SBarry Smith       }
389d32f9abdSBarry Smith 
3900298fd71SBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&fieldsplit_default,NULL);CHKERRQ(ierr);
3916ce1633cSBarry Smith       if (stokes) {
3926ce1633cSBarry Smith         IS       zerodiags,rest;
3936ce1633cSBarry Smith         PetscInt nmin,nmax;
3946ce1633cSBarry Smith 
3956ce1633cSBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
3966ce1633cSBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
3976ce1633cSBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
3987287d2a3SDmitry Karpeev         if (jac->reset) {
3997287d2a3SDmitry Karpeev           jac->head->is       = rest;
4007287d2a3SDmitry Karpeev           jac->head->next->is = zerodiags;
4012fa5cd67SKarl Rupp         } else {
4026ce1633cSBarry Smith           ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
4036ce1633cSBarry Smith           ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr);
4047287d2a3SDmitry Karpeev         }
405fcfd50ebSBarry Smith         ierr = ISDestroy(&zerodiags);CHKERRQ(ierr);
406fcfd50ebSBarry Smith         ierr = ISDestroy(&rest);CHKERRQ(ierr);
4076ce1633cSBarry Smith       } else {
408*ce94432eSBarry Smith         if (jac->reset) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used");
4097287d2a3SDmitry Karpeev         if (!fieldsplit_default) {
410d32f9abdSBarry Smith           /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
411d32f9abdSBarry Smith            then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
4126c924f48SJed Brown           ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
4136c924f48SJed Brown           if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
414d32f9abdSBarry Smith         }
4157287d2a3SDmitry Karpeev         if (fieldsplit_default || !jac->splitdefined) {
416d32f9abdSBarry Smith           ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
417db4c96c1SJed Brown           for (i=0; i<jac->bs; i++) {
4186c924f48SJed Brown             char splitname[8];
4198caf3d72SBarry Smith             ierr = PetscSNPrintf(splitname,sizeof(splitname),"%D",i);CHKERRQ(ierr);
4205d4c12cdSJungho Lee             ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr);
42179416396SBarry Smith           }
4225d4c12cdSJungho Lee           jac->defaultsplit = PETSC_TRUE;
423521d7252SBarry Smith         }
42466ffff09SJed Brown       }
4256ce1633cSBarry Smith     }
426edf189efSBarry Smith   } else if (jac->nsplits == 1) {
427edf189efSBarry Smith     if (ilink->is) {
428edf189efSBarry Smith       IS       is2;
429edf189efSBarry Smith       PetscInt nmin,nmax;
430edf189efSBarry Smith 
431edf189efSBarry Smith       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
432edf189efSBarry Smith       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
433db4c96c1SJed Brown       ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr);
434fcfd50ebSBarry Smith       ierr = ISDestroy(&is2);CHKERRQ(ierr);
4357b62db95SJungho Lee     } else SETERRQ(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
436edf189efSBarry Smith   }
437d0af7cd3SBarry Smith 
438d0af7cd3SBarry Smith 
4397b62db95SJungho Lee   if (jac->nsplits < 2) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits);
44069a612a9SBarry Smith   PetscFunctionReturn(0);
44169a612a9SBarry Smith }
44269a612a9SBarry Smith 
443514bf10dSMatthew G Knepley PetscErrorCode PetscOptionsFindPairPrefix_Private(const char pre[], const char name[], char *value[], PetscBool *flg);
444514bf10dSMatthew G Knepley 
44569a612a9SBarry Smith #undef __FUNCT__
44669a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
44769a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
44869a612a9SBarry Smith {
44969a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
45069a612a9SBarry Smith   PetscErrorCode    ierr;
4515a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
4522c9966d7SBarry Smith   PetscInt          i,nsplit;
45369a612a9SBarry Smith   MatStructure      flag = pc->flag;
4545f4ab4e1SJungho Lee   PetscBool         sorted, sorted_col;
45569a612a9SBarry Smith 
45669a612a9SBarry Smith   PetscFunctionBegin;
45769a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
45897bbdb24SBarry Smith   nsplit = jac->nsplits;
4595a9f2f41SSatish Balay   ilink  = jac->head;
46097bbdb24SBarry Smith 
46197bbdb24SBarry Smith   /* get the matrices for each split */
462704ba839SBarry Smith   if (!jac->issetup) {
4631b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
46497bbdb24SBarry Smith 
465704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
466704ba839SBarry Smith 
4675d4c12cdSJungho Lee     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
4682c9966d7SBarry Smith     if (jac->defaultsplit || !ilink->is) {
4692c9966d7SBarry Smith       if (jac->bs <= 0) jac->bs = nsplit;
4702c9966d7SBarry Smith     }
47151f519a2SBarry Smith     bs     = jac->bs;
47297bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
4731b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
4741b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
4751b9fc7fcSBarry Smith       if (jac->defaultsplit) {
476*ce94432eSBarry Smith         ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
4775f4ab4e1SJungho Lee         ierr = ISDuplicate(ilink->is,&ilink->is_col);CHKERRQ(ierr);
478704ba839SBarry Smith       } else if (!ilink->is) {
479ccb205f8SBarry Smith         if (ilink->nfields > 1) {
4805f4ab4e1SJungho Lee           PetscInt *ii,*jj,j,k,nfields = ilink->nfields,*fields = ilink->fields,*fields_col = ilink->fields_col;
4815a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
4825f4ab4e1SJungho Lee           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&jj);CHKERRQ(ierr);
4831b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
4841b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
4851b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
4865f4ab4e1SJungho Lee               jj[nfields*j + k] = rstart + bs*j + fields_col[k];
48797bbdb24SBarry Smith             }
48897bbdb24SBarry Smith           }
489*ce94432eSBarry Smith           ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr);
490*ce94432eSBarry Smith           ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),nslots*nfields,jj,PETSC_OWN_POINTER,&ilink->is_col);CHKERRQ(ierr);
491ccb205f8SBarry Smith         } else {
492*ce94432eSBarry Smith           ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
493*ce94432eSBarry Smith           ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),nslots,rstart+ilink->fields_col[0],bs,&ilink->is_col);CHKERRQ(ierr);
494ccb205f8SBarry Smith        }
4953e197d65SBarry Smith       }
496704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
4975f4ab4e1SJungho Lee       if (ilink->is_col) { ierr = ISSorted(ilink->is_col,&sorted_col);CHKERRQ(ierr); }
4985f4ab4e1SJungho Lee       if (!sorted || !sorted_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
499704ba839SBarry Smith       ilink = ilink->next;
5001b9fc7fcSBarry Smith     }
5011b9fc7fcSBarry Smith   }
5021b9fc7fcSBarry Smith 
503704ba839SBarry Smith   ilink = jac->head;
50497bbdb24SBarry Smith   if (!jac->pmat) {
505bdddcaaaSMatthew G Knepley     Vec xtmp;
506bdddcaaaSMatthew G Knepley 
5070298fd71SBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,NULL);CHKERRQ(ierr);
508cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
509bdddcaaaSMatthew G Knepley     ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
510cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
511bdddcaaaSMatthew G Knepley       MatNullSpace sp;
512bdddcaaaSMatthew G Knepley 
513a3df900dSMatthew G Knepley       /* Check for preconditioning matrix attached to IS */
514a3df900dSMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &jac->pmat[i]);CHKERRQ(ierr);
515a3df900dSMatthew G Knepley       if (jac->pmat[i]) {
516a3df900dSMatthew G Knepley         ierr = PetscObjectReference((PetscObject) jac->pmat[i]);CHKERRQ(ierr);
517a3df900dSMatthew G Knepley         if (jac->type == PC_COMPOSITE_SCHUR) {
518a3df900dSMatthew G Knepley           jac->schur_user = jac->pmat[i];
5192fa5cd67SKarl Rupp 
520a3df900dSMatthew G Knepley           ierr = PetscObjectReference((PetscObject) jac->schur_user);CHKERRQ(ierr);
521a3df900dSMatthew G Knepley         }
522a3df900dSMatthew G Knepley       } else {
5235f4ab4e1SJungho Lee         ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
524a3df900dSMatthew G Knepley       }
525bdddcaaaSMatthew G Knepley       /* create work vectors for each split */
526bdddcaaaSMatthew G Knepley       ierr     = MatGetVecs(jac->pmat[i],&jac->x[i],&jac->y[i]);CHKERRQ(ierr);
5270298fd71SBarry Smith       ilink->x = jac->x[i]; ilink->y = jac->y[i]; ilink->z = NULL;
528bdddcaaaSMatthew G Knepley       /* compute scatter contexts needed by multiplicative versions and non-default splits */
5290298fd71SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],NULL,&ilink->sctx);CHKERRQ(ierr);
530ed1f0337SMatthew G Knepley       /* Check for null space attached to IS */
531bafc1b83SMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "nullspace", (PetscObject*) &sp);CHKERRQ(ierr);
532bafc1b83SMatthew G Knepley       if (sp) {
533bafc1b83SMatthew G Knepley         ierr = MatSetNullSpace(jac->pmat[i], sp);CHKERRQ(ierr);
534bafc1b83SMatthew G Knepley       }
535ed1f0337SMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "nearnullspace", (PetscObject*) &sp);CHKERRQ(ierr);
536ed1f0337SMatthew G Knepley       if (sp) {
537ed1f0337SMatthew G Knepley         ierr = MatSetNearNullSpace(jac->pmat[i], sp);CHKERRQ(ierr);
538ed1f0337SMatthew G Knepley       }
539704ba839SBarry Smith       ilink = ilink->next;
540cf502942SBarry Smith     }
541bdddcaaaSMatthew G Knepley     ierr = VecDestroy(&xtmp);CHKERRQ(ierr);
54297bbdb24SBarry Smith   } else {
543cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
544a3df900dSMatthew G Knepley       Mat pmat;
545a3df900dSMatthew G Knepley 
546a3df900dSMatthew G Knepley       /* Check for preconditioning matrix attached to IS */
547a3df900dSMatthew G Knepley       ierr = PetscObjectQuery((PetscObject) ilink->is, "pmat", (PetscObject*) &pmat);CHKERRQ(ierr);
548a3df900dSMatthew G Knepley       if (!pmat) {
5495f4ab4e1SJungho Lee         ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
550a3df900dSMatthew G Knepley       }
551704ba839SBarry Smith       ilink = ilink->next;
552cf502942SBarry Smith     }
55397bbdb24SBarry Smith   }
554519d70e2SJed Brown   if (jac->realdiagonal) {
555519d70e2SJed Brown     ilink = jac->head;
556519d70e2SJed Brown     if (!jac->mat) {
557519d70e2SJed Brown       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr);
558519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
5595f4ab4e1SJungho Lee         ierr  = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
560519d70e2SJed Brown         ilink = ilink->next;
561519d70e2SJed Brown       }
562519d70e2SJed Brown     } else {
563519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
5645f4ab4e1SJungho Lee         if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is_col,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);}
565519d70e2SJed Brown         ilink = ilink->next;
566519d70e2SJed Brown       }
567519d70e2SJed Brown     }
568519d70e2SJed Brown   } else {
569519d70e2SJed Brown     jac->mat = jac->pmat;
570519d70e2SJed Brown   }
57197bbdb24SBarry Smith 
5726c8605c2SJed Brown   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
57368dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
57468dd23aaSBarry Smith     ilink = jac->head;
57568dd23aaSBarry Smith     if (!jac->Afield) {
57668dd23aaSBarry Smith       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
57768dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
5780298fd71SBarry Smith         ierr  = MatGetSubMatrix(pc->mat,ilink->is,NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
57968dd23aaSBarry Smith         ilink = ilink->next;
58068dd23aaSBarry Smith       }
58168dd23aaSBarry Smith     } else {
58268dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
5830298fd71SBarry Smith         ierr  = MatGetSubMatrix(pc->mat,ilink->is,NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
58468dd23aaSBarry Smith         ilink = ilink->next;
58568dd23aaSBarry Smith       }
58668dd23aaSBarry Smith     }
58768dd23aaSBarry Smith   }
58868dd23aaSBarry Smith 
5893b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
5903b224e63SBarry Smith     IS          ccis;
5914aa3045dSJed Brown     PetscInt    rstart,rend;
592093c86ffSJed Brown     char        lscname[256];
593093c86ffSJed Brown     PetscObject LSC_L;
594*ce94432eSBarry Smith 
595*ce94432eSBarry Smith     if (nsplit != 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
59668dd23aaSBarry Smith 
597e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
598e6cab6aaSJed Brown     ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
599e6cab6aaSJed Brown 
6003b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
6013b224e63SBarry Smith     if (jac->schur) {
6020298fd71SBarry Smith       KSP kspA = jac->head->ksp, kspInner = NULL, kspUpper = jac->kspupper;
603443836d0SMatthew G Knepley 
604fb3147dbSMatthew G Knepley       ierr  = MatSchurComplementGetKSP(jac->schur, &kspInner);CHKERRQ(ierr);
6053b224e63SBarry Smith       ilink = jac->head;
60649bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
6074aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
608fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
6093b224e63SBarry Smith       ilink = ilink->next;
61049bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
6114aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
612fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
613a3df900dSMatthew G Knepley       ierr  = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],pc->flag);CHKERRQ(ierr);
614443836d0SMatthew G Knepley       if (kspA != kspInner) {
615443836d0SMatthew G Knepley         ierr = KSPSetOperators(kspA,jac->mat[0],jac->pmat[0],pc->flag);CHKERRQ(ierr);
616443836d0SMatthew G Knepley       }
617443836d0SMatthew G Knepley       if (kspUpper != kspA) {
618443836d0SMatthew G Knepley         ierr = KSPSetOperators(kspUpper,jac->mat[0],jac->pmat[0],pc->flag);CHKERRQ(ierr);
619443836d0SMatthew G Knepley       }
620084e4875SJed Brown       ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr);
6213b224e63SBarry Smith     } else {
6221cee3971SBarry Smith       KSP          ksp;
623bafc1b83SMatthew G Knepley       const char   *Dprefix;
624bafc1b83SMatthew G Knepley       char         schurprefix[256];
625514bf10dSMatthew G Knepley       char         schurtestoption[256];
626bdddcaaaSMatthew G Knepley       MatNullSpace sp;
627514bf10dSMatthew G Knepley       PetscBool    flg;
6283b224e63SBarry Smith 
629a04f6461SBarry Smith       /* extract the A01 and A10 matrices */
6303b224e63SBarry Smith       ilink = jac->head;
63149bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
6324aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
633fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
6343b224e63SBarry Smith       ilink = ilink->next;
63549bb4cd7SJungho Lee       ierr  = ISComplement(ilink->is_col,rstart,rend,&ccis);CHKERRQ(ierr);
6364aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
637fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
63820252d06SBarry Smith 
63920252d06SBarry Smith       /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] to define Schur complement */
64020252d06SBarry Smith       ierr = MatCreate(((PetscObject)jac->mat[0])->comm,&jac->schur);CHKERRQ(ierr);
64120252d06SBarry Smith       ierr = MatSetType(jac->schur,MATSCHURCOMPLEMENT);CHKERRQ(ierr);
64220252d06SBarry Smith       ierr = MatSchurComplementGetKSP(jac->schur, &ksp);CHKERRQ(ierr);
64320252d06SBarry Smith       ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_inner_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr);
64420252d06SBarry Smith       /* Indent this deeper to emphasize the "inner" nature of this solver. */
64520252d06SBarry Smith       ierr = PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject) pc, 2);CHKERRQ(ierr);
64620252d06SBarry Smith       ierr = KSPSetOptionsPrefix(ksp, schurprefix);CHKERRQ(ierr);
64720252d06SBarry Smith       ierr = MatSchurComplementSet(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1]);CHKERRQ(ierr);
64820252d06SBarry Smith 
649bdddcaaaSMatthew G Knepley       ierr = MatGetNullSpace(jac->pmat[1], &sp);CHKERRQ(ierr);
65020252d06SBarry Smith       if (sp) {
65120252d06SBarry Smith         ierr = MatSetNullSpace(jac->schur, sp);CHKERRQ(ierr);
65220252d06SBarry Smith       }
65320252d06SBarry Smith 
65420252d06SBarry Smith       ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_inner_", ilink->splitname);CHKERRQ(ierr);
6550298fd71SBarry Smith       ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);CHKERRQ(ierr);
656514bf10dSMatthew G Knepley       if (flg) {
657514bf10dSMatthew G Knepley         DM dmInner;
658514bf10dSMatthew G Knepley 
659514bf10dSMatthew G Knepley         /* Set DM for new solver */
660bafc1b83SMatthew G Knepley         ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr);
661bafc1b83SMatthew G Knepley         ierr = KSPSetDM(ksp, dmInner);CHKERRQ(ierr);
6627287d2a3SDmitry Karpeev         ierr = KSPSetDMActive(ksp, PETSC_FALSE);CHKERRQ(ierr);
663443836d0SMatthew G Knepley         ierr = KSPSetOperators(jac->head->ksp,jac->mat[0],jac->pmat[0],flag);CHKERRQ(ierr);
664443836d0SMatthew G Knepley         ierr = KSPSetFromOptions(jac->head->ksp);CHKERRQ(ierr);
665514bf10dSMatthew G Knepley       } else {
666514bf10dSMatthew G Knepley         ierr = MatSchurComplementSetKSP(jac->schur, jac->head->ksp);CHKERRQ(ierr);
667514bf10dSMatthew G Knepley         ierr = PetscObjectReference((PetscObject) jac->head->ksp);CHKERRQ(ierr);
668bafc1b83SMatthew G Knepley       }
66920b26d62SBarry Smith       ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
6703b224e63SBarry Smith 
671443836d0SMatthew G Knepley       ierr = PetscSNPrintf(schurtestoption, sizeof(schurtestoption), "-fieldsplit_%s_upper_", ilink->splitname);CHKERRQ(ierr);
6720298fd71SBarry Smith       ierr = PetscOptionsFindPairPrefix_Private(((PetscObject)pc)->prefix, schurtestoption, NULL, &flg);CHKERRQ(ierr);
673443836d0SMatthew G Knepley       if (flg) {
674443836d0SMatthew G Knepley         DM dmInner;
675443836d0SMatthew G Knepley 
676443836d0SMatthew G Knepley         ierr = PetscSNPrintf(schurprefix, sizeof(schurprefix), "%sfieldsplit_%s_upper_", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ilink->splitname);CHKERRQ(ierr);
677443836d0SMatthew G Knepley         ierr = KSPCreate(((PetscObject) pc)->comm, &jac->kspupper);CHKERRQ(ierr);
678443836d0SMatthew G Knepley         ierr = KSPSetOptionsPrefix(jac->kspupper, schurprefix);CHKERRQ(ierr);
679443836d0SMatthew G Knepley         ierr = KSPGetDM(jac->head->ksp, &dmInner);CHKERRQ(ierr);
680443836d0SMatthew G Knepley         ierr = KSPSetDM(jac->kspupper, dmInner);CHKERRQ(ierr);
681443836d0SMatthew G Knepley         ierr = KSPSetDMActive(jac->kspupper, PETSC_FALSE);CHKERRQ(ierr);
682443836d0SMatthew G Knepley         ierr = KSPSetFromOptions(jac->kspupper);CHKERRQ(ierr);
683443836d0SMatthew G Knepley         ierr = KSPSetOperators(jac->kspupper,jac->mat[0],jac->pmat[0],flag);CHKERRQ(ierr);
684443836d0SMatthew G Knepley         ierr = VecDuplicate(jac->head->x, &jac->head->z);CHKERRQ(ierr);
685443836d0SMatthew G Knepley       } else {
686443836d0SMatthew G Knepley         jac->kspupper = jac->head->ksp;
687443836d0SMatthew G Knepley         ierr          = PetscObjectReference((PetscObject) jac->head->ksp);CHKERRQ(ierr);
688443836d0SMatthew G Knepley       }
689443836d0SMatthew G Knepley 
690*ce94432eSBarry Smith       ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&jac->kspschur);CHKERRQ(ierr);
6919005cf84SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr);
69220252d06SBarry Smith       ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
693084e4875SJed Brown       ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
694084e4875SJed Brown       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
6957233a360SDmitry Karpeev         PC pcschur;
6967233a360SDmitry Karpeev         ierr = KSPGetPC(jac->kspschur,&pcschur);CHKERRQ(ierr);
6977233a360SDmitry Karpeev         ierr = PCSetType(pcschur,PCNONE);CHKERRQ(ierr);
698084e4875SJed Brown         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
699e69d4d44SBarry Smith       }
7007287d2a3SDmitry Karpeev       ierr = KSPGetOptionsPrefix(jac->head->next->ksp, &Dprefix);CHKERRQ(ierr);
7017287d2a3SDmitry Karpeev       ierr = KSPSetOptionsPrefix(jac->kspschur,         Dprefix);CHKERRQ(ierr);
7023b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
70320b26d62SBarry Smith       /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */
70420b26d62SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
7053b224e63SBarry Smith     }
706093c86ffSJed Brown 
707093c86ffSJed Brown     /* HACK: special support to forward L and Lp matrices that might be used by PCLSC */
7088caf3d72SBarry Smith     ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_L",ilink->splitname);CHKERRQ(ierr);
709093c86ffSJed Brown     ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
710093c86ffSJed Brown     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
711093c86ffSJed Brown     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_L",(PetscObject)LSC_L);CHKERRQ(ierr);}
7128caf3d72SBarry Smith     ierr = PetscSNPrintf(lscname,sizeof(lscname),"%s_LSC_Lp",ilink->splitname);CHKERRQ(ierr);
713093c86ffSJed Brown     ierr = PetscObjectQuery((PetscObject)pc->pmat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);
714093c86ffSJed Brown     if (!LSC_L) {ierr = PetscObjectQuery((PetscObject)pc->mat,lscname,(PetscObject*)&LSC_L);CHKERRQ(ierr);}
715093c86ffSJed Brown     if (LSC_L) {ierr = PetscObjectCompose((PetscObject)jac->schur,"LSC_Lp",(PetscObject)LSC_L);CHKERRQ(ierr);}
7163b224e63SBarry Smith   } else {
71768bd789dSDmitry Karpeev     /* set up the individual splits' PCs */
71897bbdb24SBarry Smith     i     = 0;
7195a9f2f41SSatish Balay     ilink = jac->head;
7205a9f2f41SSatish Balay     while (ilink) {
721519d70e2SJed Brown       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr);
7223b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
723c1570756SJed Brown       if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);}
72497bbdb24SBarry Smith       i++;
7255a9f2f41SSatish Balay       ilink = ilink->next;
7260971522cSBarry Smith     }
7273b224e63SBarry Smith   }
7283b224e63SBarry Smith 
729c1570756SJed Brown   jac->suboptionsset = PETSC_TRUE;
7300971522cSBarry Smith   PetscFunctionReturn(0);
7310971522cSBarry Smith }
7320971522cSBarry Smith 
7335a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
734ca9f406cSSatish Balay   (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
735ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
7365a9f2f41SSatish Balay    KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
737ca9f406cSSatish Balay    VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
738ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
73979416396SBarry Smith 
7400971522cSBarry Smith #undef __FUNCT__
7413b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur"
7423b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
7433b224e63SBarry Smith {
7443b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
7453b224e63SBarry Smith   PetscErrorCode    ierr;
7463b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
747443836d0SMatthew G Knepley   KSP               kspA   = ilinkA->ksp, kspLower = kspA, kspUpper = jac->kspupper;
7483b224e63SBarry Smith 
7493b224e63SBarry Smith   PetscFunctionBegin;
750c5d2311dSJed Brown   switch (jac->schurfactorization) {
751c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_DIAG:
752a04f6461SBarry Smith     /* [A00 0; 0 -S], positive definite, suitable for MINRES */
753c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
754c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
755c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
756443836d0SMatthew G Knepley     ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
757c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
758c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
759c5d2311dSJed Brown     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
760c5d2311dSJed Brown     ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr);
761c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
762c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
763c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
764c5d2311dSJed Brown     break;
765c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_LOWER:
766a04f6461SBarry Smith     /* [A00 0; A10 S], suitable for left preconditioning */
767c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
768c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
769443836d0SMatthew G Knepley     ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
770c5d2311dSJed Brown     ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
771c5d2311dSJed Brown     ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr);
772c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
773c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
774c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
775c5d2311dSJed Brown     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
776c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
777c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
778c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
779c5d2311dSJed Brown     break;
780c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_UPPER:
781a04f6461SBarry Smith     /* [A00 A01; 0 S], suitable for right preconditioning */
782c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
783c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
784c5d2311dSJed Brown     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
785c5d2311dSJed Brown     ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
786c5d2311dSJed Brown     ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr);
787c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
788c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
789c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
790443836d0SMatthew G Knepley     ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
791c5d2311dSJed Brown     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
792c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
793c5d2311dSJed Brown     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
794c5d2311dSJed Brown     break;
795c9c6ffaaSJed Brown   case PC_FIELDSPLIT_SCHUR_FACT_FULL:
796a04f6461SBarry 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 */
7973b224e63SBarry Smith     ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7983b224e63SBarry Smith     ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
799443836d0SMatthew G Knepley     ierr = KSPSolve(kspLower,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
8003b224e63SBarry Smith     ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
8013b224e63SBarry Smith     ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
8023b224e63SBarry Smith     ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8033b224e63SBarry Smith     ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8043b224e63SBarry Smith 
8053b224e63SBarry Smith     ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
8063b224e63SBarry Smith     ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8073b224e63SBarry Smith     ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8083b224e63SBarry Smith 
809443836d0SMatthew G Knepley     if (kspUpper == kspA) {
8103b224e63SBarry Smith       ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
8113b224e63SBarry Smith       ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
812443836d0SMatthew G Knepley       ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
813443836d0SMatthew G Knepley     } else {
814443836d0SMatthew G Knepley       ierr = KSPSolve(kspA,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
815443836d0SMatthew G Knepley       ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
816443836d0SMatthew G Knepley       ierr = KSPSolve(kspUpper,ilinkA->x,ilinkA->z);CHKERRQ(ierr);
817443836d0SMatthew G Knepley       ierr = VecAXPY(ilinkA->y,-1.0,ilinkA->z);CHKERRQ(ierr);
818443836d0SMatthew G Knepley     }
8193b224e63SBarry Smith     ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8203b224e63SBarry Smith     ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
821c5d2311dSJed Brown   }
8223b224e63SBarry Smith   PetscFunctionReturn(0);
8233b224e63SBarry Smith }
8243b224e63SBarry Smith 
8253b224e63SBarry Smith #undef __FUNCT__
8260971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
8270971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
8280971522cSBarry Smith {
8290971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
8300971522cSBarry Smith   PetscErrorCode    ierr;
8315a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
832939b8a20SBarry Smith   PetscInt          cnt,bs;
8330971522cSBarry Smith 
8340971522cSBarry Smith   PetscFunctionBegin;
8355847835fSJed Brown   if (jac->bs > 0) {
836077aedafSJed Brown     ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
837077aedafSJed Brown     ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
8385847835fSJed Brown   }
83951f519a2SBarry Smith   CHKMEMQ;
84079416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
8411b9fc7fcSBarry Smith     if (jac->defaultsplit) {
842939b8a20SBarry Smith       ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
843*ce94432eSBarry 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);
844939b8a20SBarry Smith       ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr);
845*ce94432eSBarry 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);
8460971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
8475a9f2f41SSatish Balay       while (ilink) {
8485a9f2f41SSatish Balay         ierr  = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
8495a9f2f41SSatish Balay         ilink = ilink->next;
8500971522cSBarry Smith       }
8510971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
8521b9fc7fcSBarry Smith     } else {
853efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
8545a9f2f41SSatish Balay       while (ilink) {
8555a9f2f41SSatish Balay         ierr  = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
8565a9f2f41SSatish Balay         ilink = ilink->next;
8571b9fc7fcSBarry Smith       }
8581b9fc7fcSBarry Smith     }
85916913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
86079416396SBarry Smith     if (!jac->w1) {
86179416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
86279416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
86379416396SBarry Smith     }
864efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
8655a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
8663e197d65SBarry Smith     cnt  = 1;
8675a9f2f41SSatish Balay     while (ilink->next) {
8685a9f2f41SSatish Balay       ilink = ilink->next;
8693e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
8703e197d65SBarry Smith       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
8713e197d65SBarry Smith       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
8723e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8733e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8743e197d65SBarry Smith       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
8753e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8763e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8773e197d65SBarry Smith     }
87851f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
87911755939SBarry Smith       cnt -= 2;
88051f519a2SBarry Smith       while (ilink->previous) {
88151f519a2SBarry Smith         ilink = ilink->previous;
88211755939SBarry Smith         /* compute the residual only over the part of the vector needed */
88311755939SBarry Smith         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
88411755939SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
88511755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
88611755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
88711755939SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
88811755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
88911755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
89051f519a2SBarry Smith       }
89111755939SBarry Smith     }
892*ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
89351f519a2SBarry Smith   CHKMEMQ;
8940971522cSBarry Smith   PetscFunctionReturn(0);
8950971522cSBarry Smith }
8960971522cSBarry Smith 
897421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
898ca9f406cSSatish Balay   (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
899ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
900421e10b8SBarry Smith    KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
901ca9f406cSSatish Balay    VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
902ca9f406cSSatish Balay    VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
903421e10b8SBarry Smith 
904421e10b8SBarry Smith #undef __FUNCT__
9058c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit"
906421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
907421e10b8SBarry Smith {
908421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
909421e10b8SBarry Smith   PetscErrorCode    ierr;
910421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
911939b8a20SBarry Smith   PetscInt          bs;
912421e10b8SBarry Smith 
913421e10b8SBarry Smith   PetscFunctionBegin;
914421e10b8SBarry Smith   CHKMEMQ;
915421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
916421e10b8SBarry Smith     if (jac->defaultsplit) {
917939b8a20SBarry Smith       ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
918*ce94432eSBarry 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);
919939b8a20SBarry Smith       ierr = VecGetBlockSize(y,&bs);CHKERRQ(ierr);
920*ce94432eSBarry 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);
921421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
922421e10b8SBarry Smith       while (ilink) {
923421e10b8SBarry Smith         ierr  = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
924421e10b8SBarry Smith         ilink = ilink->next;
925421e10b8SBarry Smith       }
926421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
927421e10b8SBarry Smith     } else {
928421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
929421e10b8SBarry Smith       while (ilink) {
930421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
931421e10b8SBarry Smith         ilink = ilink->next;
932421e10b8SBarry Smith       }
933421e10b8SBarry Smith     }
934421e10b8SBarry Smith   } else {
935421e10b8SBarry Smith     if (!jac->w1) {
936421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
937421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
938421e10b8SBarry Smith     }
939421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
940421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
941421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
942421e10b8SBarry Smith       while (ilink->next) {
943421e10b8SBarry Smith         ilink = ilink->next;
9449989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
945421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
946421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
947421e10b8SBarry Smith       }
948421e10b8SBarry Smith       while (ilink->previous) {
949421e10b8SBarry Smith         ilink = ilink->previous;
9509989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
951421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
952421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
953421e10b8SBarry Smith       }
954421e10b8SBarry Smith     } else {
955421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
956421e10b8SBarry Smith         ilink = ilink->next;
957421e10b8SBarry Smith       }
958421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
959421e10b8SBarry Smith       while (ilink->previous) {
960421e10b8SBarry Smith         ilink = ilink->previous;
9619989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
962421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
963421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
964421e10b8SBarry Smith       }
965421e10b8SBarry Smith     }
966421e10b8SBarry Smith   }
967421e10b8SBarry Smith   CHKMEMQ;
968421e10b8SBarry Smith   PetscFunctionReturn(0);
969421e10b8SBarry Smith }
970421e10b8SBarry Smith 
9710971522cSBarry Smith #undef __FUNCT__
972574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit"
973574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc)
9740971522cSBarry Smith {
9750971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
9760971522cSBarry Smith   PetscErrorCode    ierr;
9775a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
9780971522cSBarry Smith 
9790971522cSBarry Smith   PetscFunctionBegin;
9805a9f2f41SSatish Balay   while (ilink) {
981574deadeSBarry Smith     ierr  = KSPReset(ilink->ksp);CHKERRQ(ierr);
982fcfd50ebSBarry Smith     ierr  = VecDestroy(&ilink->x);CHKERRQ(ierr);
983fcfd50ebSBarry Smith     ierr  = VecDestroy(&ilink->y);CHKERRQ(ierr);
984443836d0SMatthew G Knepley     ierr  = VecDestroy(&ilink->z);CHKERRQ(ierr);
985fcfd50ebSBarry Smith     ierr  = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr);
986fcfd50ebSBarry Smith     ierr  = ISDestroy(&ilink->is);CHKERRQ(ierr);
987b5787286SJed Brown     ierr  = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
9885a9f2f41SSatish Balay     next  = ilink->next;
9895a9f2f41SSatish Balay     ilink = next;
9900971522cSBarry Smith   }
99105b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
992574deadeSBarry Smith   if (jac->mat && jac->mat != jac->pmat) {
993574deadeSBarry Smith     ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);
994574deadeSBarry Smith   } else if (jac->mat) {
9950298fd71SBarry Smith     jac->mat = NULL;
996574deadeSBarry Smith   }
99797bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
99868dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
9996bf464f9SBarry Smith   ierr       = VecDestroy(&jac->w1);CHKERRQ(ierr);
10006bf464f9SBarry Smith   ierr       = VecDestroy(&jac->w2);CHKERRQ(ierr);
10016bf464f9SBarry Smith   ierr       = MatDestroy(&jac->schur);CHKERRQ(ierr);
10026bf464f9SBarry Smith   ierr       = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
10036bf464f9SBarry Smith   ierr       = KSPDestroy(&jac->kspschur);CHKERRQ(ierr);
1004d78dad28SBarry Smith   ierr       = KSPDestroy(&jac->kspupper);CHKERRQ(ierr);
10056bf464f9SBarry Smith   ierr       = MatDestroy(&jac->B);CHKERRQ(ierr);
10066bf464f9SBarry Smith   ierr       = MatDestroy(&jac->C);CHKERRQ(ierr);
100763ec74ffSBarry Smith   jac->reset = PETSC_TRUE;
1008574deadeSBarry Smith   PetscFunctionReturn(0);
1009574deadeSBarry Smith }
1010574deadeSBarry Smith 
1011574deadeSBarry Smith #undef __FUNCT__
1012574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
1013574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
1014574deadeSBarry Smith {
1015574deadeSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1016574deadeSBarry Smith   PetscErrorCode    ierr;
1017574deadeSBarry Smith   PC_FieldSplitLink ilink = jac->head,next;
1018574deadeSBarry Smith 
1019574deadeSBarry Smith   PetscFunctionBegin;
1020574deadeSBarry Smith   ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr);
1021574deadeSBarry Smith   while (ilink) {
10226bf464f9SBarry Smith     ierr  = KSPDestroy(&ilink->ksp);CHKERRQ(ierr);
1023574deadeSBarry Smith     next  = ilink->next;
1024574deadeSBarry Smith     ierr  = PetscFree(ilink->splitname);CHKERRQ(ierr);
1025574deadeSBarry Smith     ierr  = PetscFree(ilink->fields);CHKERRQ(ierr);
10265d4c12cdSJungho Lee     ierr  = PetscFree(ilink->fields_col);CHKERRQ(ierr);
1027574deadeSBarry Smith     ierr  = PetscFree(ilink);CHKERRQ(ierr);
1028574deadeSBarry Smith     ilink = next;
1029574deadeSBarry Smith   }
1030574deadeSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
1031c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
10320298fd71SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","",NULL);CHKERRQ(ierr);
10330298fd71SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","",NULL);CHKERRQ(ierr);
10340298fd71SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","",NULL);CHKERRQ(ierr);
10350298fd71SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","",NULL);CHKERRQ(ierr);
10360298fd71SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","",NULL);CHKERRQ(ierr);
10370298fd71SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",NULL);CHKERRQ(ierr);
10380298fd71SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",NULL);CHKERRQ(ierr);
10390971522cSBarry Smith   PetscFunctionReturn(0);
10400971522cSBarry Smith }
10410971522cSBarry Smith 
10420971522cSBarry Smith #undef __FUNCT__
10430971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
10440971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
10450971522cSBarry Smith {
10461b9fc7fcSBarry Smith   PetscErrorCode  ierr;
10476c924f48SJed Brown   PetscInt        bs;
1048bc59fbc5SBarry Smith   PetscBool       flg,stokes = PETSC_FALSE;
10499dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
10503b224e63SBarry Smith   PCCompositeType ctype;
10511b9fc7fcSBarry Smith 
10520971522cSBarry Smith   PetscFunctionBegin;
10531b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
10540298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,NULL);CHKERRQ(ierr);
105551f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
105651f519a2SBarry Smith   if (flg) {
105751f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
105851f519a2SBarry Smith   }
1059704ba839SBarry Smith 
10600298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,NULL);CHKERRQ(ierr);
1061c0adfefeSBarry Smith   if (stokes) {
1062c0adfefeSBarry Smith     ierr          = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
1063c0adfefeSBarry Smith     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
1064c0adfefeSBarry Smith   }
1065c0adfefeSBarry Smith 
10663b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
10673b224e63SBarry Smith   if (flg) {
10683b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
10693b224e63SBarry Smith   }
1070c30613efSMatthew Knepley   /* Only setup fields once */
1071c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
1072d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
1073d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
10746c924f48SJed Brown     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
10756c924f48SJed Brown     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
1076d32f9abdSBarry Smith   }
1077c5d2311dSJed Brown   if (jac->type == PC_COMPOSITE_SCHUR) {
1078c9c6ffaaSJed Brown     ierr = PetscOptionsGetEnum(((PetscObject)pc)->prefix,"-pc_fieldsplit_schur_factorization_type",PCFieldSplitSchurFactTypes,(PetscEnum*)&jac->schurfactorization,&flg);CHKERRQ(ierr);
1079c9c6ffaaSJed Brown     if (flg) {ierr = PetscInfo(pc,"Deprecated use of -pc_fieldsplit_schur_factorization_type\n");CHKERRQ(ierr);}
10800298fd71SBarry 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);
10810298fd71SBarry 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);
1082c5d2311dSJed Brown   }
10831b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
10840971522cSBarry Smith   PetscFunctionReturn(0);
10850971522cSBarry Smith }
10860971522cSBarry Smith 
10870971522cSBarry Smith /*------------------------------------------------------------------------------------*/
10880971522cSBarry Smith 
10890971522cSBarry Smith EXTERN_C_BEGIN
10900971522cSBarry Smith #undef __FUNCT__
10910971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
10925d4c12cdSJungho Lee PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
10930971522cSBarry Smith {
109497bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
10950971522cSBarry Smith   PetscErrorCode    ierr;
10965a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
109769a612a9SBarry Smith   char              prefix[128];
10985d4c12cdSJungho Lee   PetscInt          i;
10990971522cSBarry Smith 
11000971522cSBarry Smith   PetscFunctionBegin;
11016c924f48SJed Brown   if (jac->splitdefined) {
11026c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
11036c924f48SJed Brown     PetscFunctionReturn(0);
11046c924f48SJed Brown   }
110551f519a2SBarry Smith   for (i=0; i<n; i++) {
1106e32f2f54SBarry 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);
1107e32f2f54SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
110851f519a2SBarry Smith   }
1109704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
1110a04f6461SBarry Smith   if (splitname) {
1111db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1112a04f6461SBarry Smith   } else {
1113a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
1114a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
1115a04f6461SBarry Smith   }
1116704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
11175a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
11185d4c12cdSJungho Lee   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields_col);CHKERRQ(ierr);
11195d4c12cdSJungho Lee   ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr);
11202fa5cd67SKarl Rupp 
11215a9f2f41SSatish Balay   ilink->nfields = n;
11220298fd71SBarry Smith   ilink->next    = NULL;
1123*ce94432eSBarry Smith   ierr           = KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);CHKERRQ(ierr);
112420252d06SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
11255a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
11269005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
112769a612a9SBarry Smith 
11288caf3d72SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);CHKERRQ(ierr);
11295a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
11300971522cSBarry Smith 
11310971522cSBarry Smith   if (!next) {
11325a9f2f41SSatish Balay     jac->head       = ilink;
11330298fd71SBarry Smith     ilink->previous = NULL;
11340971522cSBarry Smith   } else {
11350971522cSBarry Smith     while (next->next) {
11360971522cSBarry Smith       next = next->next;
11370971522cSBarry Smith     }
11385a9f2f41SSatish Balay     next->next      = ilink;
113951f519a2SBarry Smith     ilink->previous = next;
11400971522cSBarry Smith   }
11410971522cSBarry Smith   jac->nsplits++;
11420971522cSBarry Smith   PetscFunctionReturn(0);
11430971522cSBarry Smith }
11440971522cSBarry Smith EXTERN_C_END
11450971522cSBarry Smith 
1146e69d4d44SBarry Smith EXTERN_C_BEGIN
1147e69d4d44SBarry Smith #undef __FUNCT__
1148e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
11497087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
1150e69d4d44SBarry Smith {
1151e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1152e69d4d44SBarry Smith   PetscErrorCode ierr;
1153e69d4d44SBarry Smith 
1154e69d4d44SBarry Smith   PetscFunctionBegin;
1155519d70e2SJed Brown   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
1156e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
11572fa5cd67SKarl Rupp 
1158e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
115913e0d083SBarry Smith   if (n) *n = jac->nsplits;
1160e69d4d44SBarry Smith   PetscFunctionReturn(0);
1161e69d4d44SBarry Smith }
1162e69d4d44SBarry Smith EXTERN_C_END
11630971522cSBarry Smith 
11640971522cSBarry Smith EXTERN_C_BEGIN
11650971522cSBarry Smith #undef __FUNCT__
116669a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
11677087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
11680971522cSBarry Smith {
11690971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
11700971522cSBarry Smith   PetscErrorCode    ierr;
11710971522cSBarry Smith   PetscInt          cnt   = 0;
11725a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
11730971522cSBarry Smith 
11740971522cSBarry Smith   PetscFunctionBegin;
11755d480477SMatthew G Knepley   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
11765a9f2f41SSatish Balay   while (ilink) {
11775a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
11785a9f2f41SSatish Balay     ilink            = ilink->next;
11790971522cSBarry Smith   }
11805d480477SMatthew 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);
118113e0d083SBarry Smith   if (n) *n = jac->nsplits;
11820971522cSBarry Smith   PetscFunctionReturn(0);
11830971522cSBarry Smith }
11840971522cSBarry Smith EXTERN_C_END
11850971522cSBarry Smith 
1186704ba839SBarry Smith EXTERN_C_BEGIN
1187704ba839SBarry Smith #undef __FUNCT__
1188704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
11897087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
1190704ba839SBarry Smith {
1191704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1192704ba839SBarry Smith   PetscErrorCode    ierr;
1193704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
1194704ba839SBarry Smith   char              prefix[128];
1195704ba839SBarry Smith 
1196704ba839SBarry Smith   PetscFunctionBegin;
11976c924f48SJed Brown   if (jac->splitdefined) {
11986c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
11996c924f48SJed Brown     PetscFunctionReturn(0);
12006c924f48SJed Brown   }
120116913363SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
1202a04f6461SBarry Smith   if (splitname) {
1203db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1204a04f6461SBarry Smith   } else {
1205b5787286SJed Brown     ierr = PetscMalloc(8*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
1206b5787286SJed Brown     ierr = PetscSNPrintf(ilink->splitname,7,"%D",jac->nsplits);CHKERRQ(ierr);
1207a04f6461SBarry Smith   }
12081ab39975SBarry Smith   ierr          = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1209b5787286SJed Brown   ierr          = ISDestroy(&ilink->is);CHKERRQ(ierr);
1210b5787286SJed Brown   ilink->is     = is;
1211b5787286SJed Brown   ierr          = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1212b5787286SJed Brown   ierr          = ISDestroy(&ilink->is_col);CHKERRQ(ierr);
1213b5787286SJed Brown   ilink->is_col = is;
12140298fd71SBarry Smith   ilink->next   = NULL;
1215*ce94432eSBarry Smith   ierr          = KSPCreate(PetscObjectComm((PetscObject)pc),&ilink->ksp);CHKERRQ(ierr);
121620252d06SBarry Smith   ierr          = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1217704ba839SBarry Smith   ierr          = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
12189005cf84SBarry Smith   ierr          = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
1219704ba839SBarry Smith 
12208caf3d72SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%sfieldsplit_%s_",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",ilink->splitname);CHKERRQ(ierr);
1221704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
1222704ba839SBarry Smith 
1223704ba839SBarry Smith   if (!next) {
1224704ba839SBarry Smith     jac->head       = ilink;
12250298fd71SBarry Smith     ilink->previous = NULL;
1226704ba839SBarry Smith   } else {
1227704ba839SBarry Smith     while (next->next) {
1228704ba839SBarry Smith       next = next->next;
1229704ba839SBarry Smith     }
1230704ba839SBarry Smith     next->next      = ilink;
1231704ba839SBarry Smith     ilink->previous = next;
1232704ba839SBarry Smith   }
1233704ba839SBarry Smith   jac->nsplits++;
1234704ba839SBarry Smith   PetscFunctionReturn(0);
1235704ba839SBarry Smith }
1236704ba839SBarry Smith EXTERN_C_END
1237704ba839SBarry Smith 
12380971522cSBarry Smith #undef __FUNCT__
12390971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
12400971522cSBarry Smith /*@
12410971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
12420971522cSBarry Smith 
1243ad4df100SBarry Smith     Logically Collective on PC
12440971522cSBarry Smith 
12450971522cSBarry Smith     Input Parameters:
12460971522cSBarry Smith +   pc  - the preconditioner context
12470298fd71SBarry Smith .   splitname - name of this split, if NULL the number of the split is used
12480971522cSBarry Smith .   n - the number of fields in this split
1249db4c96c1SJed Brown -   fields - the fields in this split
12500971522cSBarry Smith 
12510971522cSBarry Smith     Level: intermediate
12520971522cSBarry Smith 
1253d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
1254d32f9abdSBarry Smith 
12557287d2a3SDmitry Karpeev      The PCFieldSplitSetFields() is for defining fields as strided blocks. For example, if the block
1256d32f9abdSBarry 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
1257d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1258d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
1259d32f9abdSBarry Smith 
1260db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
1261db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1262db4c96c1SJed Brown 
12635d4c12cdSJungho Lee      Developer Note: This routine does not actually create the IS representing the split, that is delayed
12645d4c12cdSJungho Lee      until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be
12655d4c12cdSJungho Lee      available when this routine is called.
12665d4c12cdSJungho Lee 
1267d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
12680971522cSBarry Smith 
12690971522cSBarry Smith @*/
12705d4c12cdSJungho Lee PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
12710971522cSBarry Smith {
12724ac538c5SBarry Smith   PetscErrorCode ierr;
12730971522cSBarry Smith 
12740971522cSBarry Smith   PetscFunctionBegin;
12750700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1276db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1277*ce94432eSBarry Smith   if (n < 1) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1278db4c96c1SJed Brown   PetscValidIntPointer(fields,3);
12795d4c12cdSJungho Lee   ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt*,const PetscInt*),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr);
12800971522cSBarry Smith   PetscFunctionReturn(0);
12810971522cSBarry Smith }
12820971522cSBarry Smith 
12830971522cSBarry Smith #undef __FUNCT__
1284704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
1285704ba839SBarry Smith /*@
1286704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
1287704ba839SBarry Smith 
1288ad4df100SBarry Smith     Logically Collective on PC
1289704ba839SBarry Smith 
1290704ba839SBarry Smith     Input Parameters:
1291704ba839SBarry Smith +   pc  - the preconditioner context
12920298fd71SBarry Smith .   splitname - name of this split, if NULL the number of the split is used
1293db4c96c1SJed Brown -   is - the index set that defines the vector elements in this field
1294704ba839SBarry Smith 
1295d32f9abdSBarry Smith 
1296a6ffb8dbSJed Brown     Notes:
1297a6ffb8dbSJed Brown     Use PCFieldSplitSetFields(), for fields defined by strided types.
1298a6ffb8dbSJed Brown 
1299db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
1300db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1301d32f9abdSBarry Smith 
1302704ba839SBarry Smith     Level: intermediate
1303704ba839SBarry Smith 
1304704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
1305704ba839SBarry Smith 
1306704ba839SBarry Smith @*/
13077087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1308704ba839SBarry Smith {
13094ac538c5SBarry Smith   PetscErrorCode ierr;
1310704ba839SBarry Smith 
1311704ba839SBarry Smith   PetscFunctionBegin;
13120700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13137b62db95SJungho Lee   if (splitname) PetscValidCharPointer(splitname,2);
1314db4c96c1SJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,3);
13154ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr);
1316704ba839SBarry Smith   PetscFunctionReturn(0);
1317704ba839SBarry Smith }
1318704ba839SBarry Smith 
1319704ba839SBarry Smith #undef __FUNCT__
132057a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS"
132157a9adfeSMatthew G Knepley /*@
132257a9adfeSMatthew G Knepley     PCFieldSplitGetIS - Retrieves the elements for a field as an IS
132357a9adfeSMatthew G Knepley 
132457a9adfeSMatthew G Knepley     Logically Collective on PC
132557a9adfeSMatthew G Knepley 
132657a9adfeSMatthew G Knepley     Input Parameters:
132757a9adfeSMatthew G Knepley +   pc  - the preconditioner context
132857a9adfeSMatthew G Knepley -   splitname - name of this split
132957a9adfeSMatthew G Knepley 
133057a9adfeSMatthew G Knepley     Output Parameter:
13310298fd71SBarry Smith -   is - the index set that defines the vector elements in this field, or NULL if the field is not found
133257a9adfeSMatthew G Knepley 
133357a9adfeSMatthew G Knepley     Level: intermediate
133457a9adfeSMatthew G Knepley 
133557a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS()
133657a9adfeSMatthew G Knepley 
133757a9adfeSMatthew G Knepley @*/
133857a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
133957a9adfeSMatthew G Knepley {
134057a9adfeSMatthew G Knepley   PetscErrorCode ierr;
134157a9adfeSMatthew G Knepley 
134257a9adfeSMatthew G Knepley   PetscFunctionBegin;
134357a9adfeSMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
134457a9adfeSMatthew G Knepley   PetscValidCharPointer(splitname,2);
134557a9adfeSMatthew G Knepley   PetscValidPointer(is,3);
134657a9adfeSMatthew G Knepley   {
134757a9adfeSMatthew G Knepley     PC_FieldSplit     *jac  = (PC_FieldSplit*) pc->data;
134857a9adfeSMatthew G Knepley     PC_FieldSplitLink ilink = jac->head;
134957a9adfeSMatthew G Knepley     PetscBool         found;
135057a9adfeSMatthew G Knepley 
13510298fd71SBarry Smith     *is = NULL;
135257a9adfeSMatthew G Knepley     while (ilink) {
135357a9adfeSMatthew G Knepley       ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr);
135457a9adfeSMatthew G Knepley       if (found) {
135557a9adfeSMatthew G Knepley         *is = ilink->is;
135657a9adfeSMatthew G Knepley         break;
135757a9adfeSMatthew G Knepley       }
135857a9adfeSMatthew G Knepley       ilink = ilink->next;
135957a9adfeSMatthew G Knepley     }
136057a9adfeSMatthew G Knepley   }
136157a9adfeSMatthew G Knepley   PetscFunctionReturn(0);
136257a9adfeSMatthew G Knepley }
136357a9adfeSMatthew G Knepley 
136457a9adfeSMatthew G Knepley #undef __FUNCT__
136551f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
136651f519a2SBarry Smith /*@
136751f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
136851f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
136951f519a2SBarry Smith 
1370ad4df100SBarry Smith     Logically Collective on PC
137151f519a2SBarry Smith 
137251f519a2SBarry Smith     Input Parameters:
137351f519a2SBarry Smith +   pc  - the preconditioner context
137451f519a2SBarry Smith -   bs - the block size
137551f519a2SBarry Smith 
137651f519a2SBarry Smith     Level: intermediate
137751f519a2SBarry Smith 
137851f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
137951f519a2SBarry Smith 
138051f519a2SBarry Smith @*/
13817087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
138251f519a2SBarry Smith {
13834ac538c5SBarry Smith   PetscErrorCode ierr;
138451f519a2SBarry Smith 
138551f519a2SBarry Smith   PetscFunctionBegin;
13860700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1387c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,bs,2);
13884ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr);
138951f519a2SBarry Smith   PetscFunctionReturn(0);
139051f519a2SBarry Smith }
139151f519a2SBarry Smith 
139251f519a2SBarry Smith #undef __FUNCT__
139369a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
13940971522cSBarry Smith /*@C
139569a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
13960971522cSBarry Smith 
139769a612a9SBarry Smith    Collective on KSP
13980971522cSBarry Smith 
13990971522cSBarry Smith    Input Parameter:
14000971522cSBarry Smith .  pc - the preconditioner context
14010971522cSBarry Smith 
14020971522cSBarry Smith    Output Parameters:
140313e0d083SBarry Smith +  n - the number of splits
140469a612a9SBarry Smith -  pc - the array of KSP contexts
14050971522cSBarry Smith 
14060971522cSBarry Smith    Note:
1407d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
1408d32f9abdSBarry Smith    (not the KSP just the array that contains them).
14090971522cSBarry Smith 
141069a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
14110971522cSBarry Smith 
1412196cc216SBarry Smith    Fortran Usage: You must pass in a KSP array that is large enough to contain all the local KSPs.
14130298fd71SBarry Smith       You can call PCFieldSplitGetSubKSP(pc,n,NULL_OBJECT,ierr) to determine how large the
1414196cc216SBarry Smith       KSP array must be.
1415196cc216SBarry Smith 
1416196cc216SBarry Smith 
14170971522cSBarry Smith    Level: advanced
14180971522cSBarry Smith 
14190971522cSBarry Smith .seealso: PCFIELDSPLIT
14200971522cSBarry Smith @*/
14217087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
14220971522cSBarry Smith {
14234ac538c5SBarry Smith   PetscErrorCode ierr;
14240971522cSBarry Smith 
14250971522cSBarry Smith   PetscFunctionBegin;
14260700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
142713e0d083SBarry Smith   if (n) PetscValidIntPointer(n,2);
14284ac538c5SBarry Smith   ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
14290971522cSBarry Smith   PetscFunctionReturn(0);
14300971522cSBarry Smith }
14310971522cSBarry Smith 
1432e69d4d44SBarry Smith #undef __FUNCT__
1433e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition"
1434e69d4d44SBarry Smith /*@
1435e69d4d44SBarry Smith     PCFieldSplitSchurPrecondition -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
1436a04f6461SBarry Smith       A11 matrix. Otherwise no preconditioner is used.
1437e69d4d44SBarry Smith 
1438e69d4d44SBarry Smith     Collective on PC
1439e69d4d44SBarry Smith 
1440e69d4d44SBarry Smith     Input Parameters:
1441e69d4d44SBarry Smith +   pc  - the preconditioner context
1442e87fab1bSBarry Smith .   ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_A11 (diag) is default
14430298fd71SBarry Smith -   userpre - matrix to use for preconditioning, or NULL
1444084e4875SJed Brown 
1445e69d4d44SBarry Smith     Options Database:
1446e87fab1bSBarry Smith .     -pc_fieldsplit_schur_precondition <self,user,a11> default is a11
1447e69d4d44SBarry Smith 
1448fd1303e9SJungho Lee     Notes:
1449fd1303e9SJungho Lee $    If ptype is
1450fd1303e9SJungho Lee $        user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument
1451fd1303e9SJungho Lee $             to this function).
1452e87fab1bSBarry Smith $        a11 then the preconditioner for the Schur complement is generated by the block diagonal part of the original
1453fd1303e9SJungho Lee $             matrix associated with the Schur complement (i.e. A11)
1454fd1303e9SJungho Lee $        self the preconditioner for the Schur complement is generated from the Schur complement matrix itself:
1455fd1303e9SJungho Lee $             The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC
1456fd1303e9SJungho Lee $             preconditioner
1457fd1303e9SJungho Lee 
1458e87fab1bSBarry Smith      When solving a saddle point problem, where the A11 block is identically zero, using a11 as the ptype only makes sense
1459fd1303e9SJungho Lee     with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
1460fd1303e9SJungho Lee     -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement.
1461fd1303e9SJungho Lee 
1462fd1303e9SJungho Lee     Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in
1463fd1303e9SJungho Lee      the name since it sets a proceedure to use.
1464fd1303e9SJungho Lee 
1465e69d4d44SBarry Smith     Level: intermediate
1466e69d4d44SBarry Smith 
1467fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC
1468e69d4d44SBarry Smith 
1469e69d4d44SBarry Smith @*/
14707087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1471e69d4d44SBarry Smith {
14724ac538c5SBarry Smith   PetscErrorCode ierr;
1473e69d4d44SBarry Smith 
1474e69d4d44SBarry Smith   PetscFunctionBegin;
14750700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
14764ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr);
1477e69d4d44SBarry Smith   PetscFunctionReturn(0);
1478e69d4d44SBarry Smith }
1479e69d4d44SBarry Smith 
1480e69d4d44SBarry Smith EXTERN_C_BEGIN
1481e69d4d44SBarry Smith #undef __FUNCT__
1482e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
14837087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1484e69d4d44SBarry Smith {
1485e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1486084e4875SJed Brown   PetscErrorCode ierr;
1487e69d4d44SBarry Smith 
1488e69d4d44SBarry Smith   PetscFunctionBegin;
1489084e4875SJed Brown   jac->schurpre = ptype;
1490084e4875SJed Brown   if (pre) {
14916bf464f9SBarry Smith     ierr            = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
1492084e4875SJed Brown     jac->schur_user = pre;
1493084e4875SJed Brown     ierr            = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1494084e4875SJed Brown   }
1495e69d4d44SBarry Smith   PetscFunctionReturn(0);
1496e69d4d44SBarry Smith }
1497e69d4d44SBarry Smith EXTERN_C_END
1498e69d4d44SBarry Smith 
149930ad9308SMatthew Knepley #undef __FUNCT__
1500c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType"
1501ab1df9f5SJed Brown /*@
1502c9c6ffaaSJed Brown     PCFieldSplitSetSchurFactType -  sets which blocks of the approximate block factorization to retain
1503ab1df9f5SJed Brown 
1504ab1df9f5SJed Brown     Collective on PC
1505ab1df9f5SJed Brown 
1506ab1df9f5SJed Brown     Input Parameters:
1507ab1df9f5SJed Brown +   pc  - the preconditioner context
1508c9c6ffaaSJed Brown -   ftype - which blocks of factorization to retain, PC_FIELDSPLIT_SCHUR_FACT_FULL is default
1509ab1df9f5SJed Brown 
1510ab1df9f5SJed Brown     Options Database:
1511c9c6ffaaSJed Brown .     -pc_fieldsplit_schur_fact_type <diag,lower,upper,full> default is full
1512ab1df9f5SJed Brown 
1513ab1df9f5SJed Brown 
1514ab1df9f5SJed Brown     Level: intermediate
1515ab1df9f5SJed Brown 
1516ab1df9f5SJed Brown     Notes:
1517ab1df9f5SJed Brown     The FULL factorization is
1518ab1df9f5SJed Brown 
1519ab1df9f5SJed Brown $   (A   B)  = (1       0) (A   0) (1  Ainv*B)
1520ab1df9f5SJed Brown $   (C   D)    (C*Ainv  1) (0   S) (0     1  )
1521ab1df9f5SJed Brown 
15226be4592eSBarry 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,
15236be4592eSBarry 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).
1524ab1df9f5SJed Brown 
15256be4592eSBarry 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
15266be4592eSBarry 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
15276be4592eSBarry 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,
15286be4592eSBarry 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.
1529ab1df9f5SJed Brown 
15306be4592eSBarry 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
15316be4592eSBarry 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).
1532ab1df9f5SJed Brown 
1533ab1df9f5SJed Brown     References:
1534ab1df9f5SJed Brown     Murphy, Golub, and Wathen, A note on preconditioning indefinite linear systems, SIAM J. Sci. Comput., 21 (2000) pp. 1969-1972.
1535ab1df9f5SJed Brown 
1536ab1df9f5SJed Brown     Ipsen, A note on preconditioning nonsymmetric matrices, SIAM J. Sci. Comput., 23 (2001), pp. 1050-1051.
1537ab1df9f5SJed Brown 
1538ab1df9f5SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType
1539ab1df9f5SJed Brown @*/
1540c9c6ffaaSJed Brown PetscErrorCode  PCFieldSplitSetSchurFactType(PC pc,PCFieldSplitSchurFactType ftype)
1541ab1df9f5SJed Brown {
1542ab1df9f5SJed Brown   PetscErrorCode ierr;
1543ab1df9f5SJed Brown 
1544ab1df9f5SJed Brown   PetscFunctionBegin;
1545ab1df9f5SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1546c9c6ffaaSJed Brown   ierr = PetscTryMethod(pc,"PCFieldSplitSetSchurFactType_C",(PC,PCFieldSplitSchurFactType),(pc,ftype));CHKERRQ(ierr);
1547ab1df9f5SJed Brown   PetscFunctionReturn(0);
1548ab1df9f5SJed Brown }
1549ab1df9f5SJed Brown 
1550ab1df9f5SJed Brown #undef __FUNCT__
1551c9c6ffaaSJed Brown #define __FUNCT__ "PCFieldSplitSetSchurFactType_FieldSplit"
1552c9c6ffaaSJed Brown PETSC_EXTERN_C PetscErrorCode PCFieldSplitSetSchurFactType_FieldSplit(PC pc,PCFieldSplitSchurFactType ftype)
1553ab1df9f5SJed Brown {
1554ab1df9f5SJed Brown   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1555ab1df9f5SJed Brown 
1556ab1df9f5SJed Brown   PetscFunctionBegin;
1557ab1df9f5SJed Brown   jac->schurfactorization = ftype;
1558ab1df9f5SJed Brown   PetscFunctionReturn(0);
1559ab1df9f5SJed Brown }
1560ab1df9f5SJed Brown 
1561ab1df9f5SJed Brown #undef __FUNCT__
156230ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
156330ad9308SMatthew Knepley /*@C
15648c03b21aSDmitry Karpeev    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
156530ad9308SMatthew Knepley 
156630ad9308SMatthew Knepley    Collective on KSP
156730ad9308SMatthew Knepley 
156830ad9308SMatthew Knepley    Input Parameter:
156930ad9308SMatthew Knepley .  pc - the preconditioner context
157030ad9308SMatthew Knepley 
157130ad9308SMatthew Knepley    Output Parameters:
1572a04f6461SBarry Smith +  A00 - the (0,0) block
1573a04f6461SBarry Smith .  A01 - the (0,1) block
1574a04f6461SBarry Smith .  A10 - the (1,0) block
1575a04f6461SBarry Smith -  A11 - the (1,1) block
157630ad9308SMatthew Knepley 
157730ad9308SMatthew Knepley    Level: advanced
157830ad9308SMatthew Knepley 
157930ad9308SMatthew Knepley .seealso: PCFIELDSPLIT
158030ad9308SMatthew Knepley @*/
1581a04f6461SBarry Smith PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
158230ad9308SMatthew Knepley {
158330ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;
158430ad9308SMatthew Knepley 
158530ad9308SMatthew Knepley   PetscFunctionBegin;
15860700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1587*ce94432eSBarry Smith   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
1588a04f6461SBarry Smith   if (A00) *A00 = jac->pmat[0];
1589a04f6461SBarry Smith   if (A01) *A01 = jac->B;
1590a04f6461SBarry Smith   if (A10) *A10 = jac->C;
1591a04f6461SBarry Smith   if (A11) *A11 = jac->pmat[1];
159230ad9308SMatthew Knepley   PetscFunctionReturn(0);
159330ad9308SMatthew Knepley }
159430ad9308SMatthew Knepley 
159579416396SBarry Smith EXTERN_C_BEGIN
159679416396SBarry Smith #undef __FUNCT__
159779416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
15987087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
159979416396SBarry Smith {
160079416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1601e69d4d44SBarry Smith   PetscErrorCode ierr;
160279416396SBarry Smith 
160379416396SBarry Smith   PetscFunctionBegin;
160479416396SBarry Smith   jac->type = type;
16053b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
16063b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
16073b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
16082fa5cd67SKarl Rupp 
1609e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1610e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1611c9c6ffaaSJed Brown     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","PCFieldSplitSetSchurFactType_FieldSplit",PCFieldSplitSetSchurFactType_FieldSplit);CHKERRQ(ierr);
1612e69d4d44SBarry Smith 
16133b224e63SBarry Smith   } else {
16143b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
16153b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
16162fa5cd67SKarl Rupp 
1617e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
16189e7d6b0aSBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr);
1619c9c6ffaaSJed Brown     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetSchurFactType_C","",0);CHKERRQ(ierr);
16203b224e63SBarry Smith   }
162179416396SBarry Smith   PetscFunctionReturn(0);
162279416396SBarry Smith }
162379416396SBarry Smith EXTERN_C_END
162479416396SBarry Smith 
162551f519a2SBarry Smith EXTERN_C_BEGIN
162651f519a2SBarry Smith #undef __FUNCT__
162751f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
16287087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
162951f519a2SBarry Smith {
163051f519a2SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
163151f519a2SBarry Smith 
163251f519a2SBarry Smith   PetscFunctionBegin;
1633*ce94432eSBarry Smith   if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
1634*ce94432eSBarry 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);
163551f519a2SBarry Smith   jac->bs = bs;
163651f519a2SBarry Smith   PetscFunctionReturn(0);
163751f519a2SBarry Smith }
163851f519a2SBarry Smith EXTERN_C_END
163951f519a2SBarry Smith 
164079416396SBarry Smith #undef __FUNCT__
164179416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
1642bc08b0f1SBarry Smith /*@
164379416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
164479416396SBarry Smith 
164579416396SBarry Smith    Collective on PC
164679416396SBarry Smith 
164779416396SBarry Smith    Input Parameter:
164879416396SBarry Smith .  pc - the preconditioner context
164981540f2fSBarry Smith .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
165079416396SBarry Smith 
165179416396SBarry Smith    Options Database Key:
1652a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
165379416396SBarry Smith 
1654b02e2d75SMatthew G Knepley    Level: Intermediate
165579416396SBarry Smith 
165679416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
165779416396SBarry Smith 
165879416396SBarry Smith .seealso: PCCompositeSetType()
165979416396SBarry Smith 
166079416396SBarry Smith @*/
16617087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
166279416396SBarry Smith {
16634ac538c5SBarry Smith   PetscErrorCode ierr;
166479416396SBarry Smith 
166579416396SBarry Smith   PetscFunctionBegin;
16660700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
16674ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr);
166879416396SBarry Smith   PetscFunctionReturn(0);
166979416396SBarry Smith }
167079416396SBarry Smith 
1671b02e2d75SMatthew G Knepley #undef __FUNCT__
1672b02e2d75SMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetType"
1673b02e2d75SMatthew G Knepley /*@
1674b02e2d75SMatthew G Knepley   PCFieldSplitGetType - Gets the type of fieldsplit preconditioner.
1675b02e2d75SMatthew G Knepley 
1676b02e2d75SMatthew G Knepley   Not collective
1677b02e2d75SMatthew G Knepley 
1678b02e2d75SMatthew G Knepley   Input Parameter:
1679b02e2d75SMatthew G Knepley . pc - the preconditioner context
1680b02e2d75SMatthew G Knepley 
1681b02e2d75SMatthew G Knepley   Output Parameter:
1682b02e2d75SMatthew G Knepley . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
1683b02e2d75SMatthew G Knepley 
1684b02e2d75SMatthew G Knepley   Level: Intermediate
1685b02e2d75SMatthew G Knepley 
1686b02e2d75SMatthew G Knepley .keywords: PC, set, type, composite preconditioner, additive, multiplicative
1687b02e2d75SMatthew G Knepley .seealso: PCCompositeSetType()
1688b02e2d75SMatthew G Knepley @*/
1689b02e2d75SMatthew G Knepley PetscErrorCode PCFieldSplitGetType(PC pc, PCCompositeType *type)
1690b02e2d75SMatthew G Knepley {
1691b02e2d75SMatthew G Knepley   PC_FieldSplit *jac = (PC_FieldSplit*) pc->data;
1692b02e2d75SMatthew G Knepley 
1693b02e2d75SMatthew G Knepley   PetscFunctionBegin;
1694b02e2d75SMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1695b02e2d75SMatthew G Knepley   PetscValidIntPointer(type,2);
1696b02e2d75SMatthew G Knepley   *type = jac->type;
1697b02e2d75SMatthew G Knepley   PetscFunctionReturn(0);
1698b02e2d75SMatthew G Knepley }
1699b02e2d75SMatthew G Knepley 
17000971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
17010971522cSBarry Smith /*MC
1702a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
1703a04f6461SBarry Smith                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.
17040971522cSBarry Smith 
1705edf189efSBarry Smith      To set options on the solvers for each block append -fieldsplit_ to all the PC
1706edf189efSBarry Smith         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
17070971522cSBarry Smith 
1708a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
170969a612a9SBarry Smith          and set the options directly on the resulting KSP object
17100971522cSBarry Smith 
17110971522cSBarry Smith    Level: intermediate
17120971522cSBarry Smith 
171379416396SBarry Smith    Options Database Keys:
171481540f2fSBarry Smith +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
171581540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
171681540f2fSBarry Smith                               been supplied explicitly by -pc_fieldsplit_%d_fields
171781540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
17180f188ba9SJed Brown .   -pc_fieldsplit_type <additive,multiplicative,symmetric_multiplicative,schur> - type of relaxation or factorization splitting
1719e87fab1bSBarry Smith .   -pc_fieldsplit_schur_precondition <self,user,a11> - default is a11
1720435f959eSBarry Smith .   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver
172179416396SBarry Smith 
17225d4c12cdSJungho Lee -    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
17235d4c12cdSJungho Lee      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
17245d4c12cdSJungho Lee 
1725c8a0d604SMatthew G Knepley    Notes:
1726c8a0d604SMatthew G Knepley     Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1727d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
1728d32f9abdSBarry Smith 
1729d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
1730d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1731d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
1732d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
1733d32f9abdSBarry Smith 
1734c8a0d604SMatthew G Knepley $     For the Schur complement preconditioner if J = ( A00 A01 )
1735c8a0d604SMatthew G Knepley $                                                    ( A10 A11 )
1736c8a0d604SMatthew G Knepley $     the preconditioner using full factorization is
1737c8a0d604SMatthew G Knepley $              ( I   -A10 ksp(A00) ) ( inv(A00)     0  ) (     I          0  )
1738c8a0d604SMatthew G Knepley $              ( 0         I       ) (   0      ksp(S) ) ( -A10 ksp(A00)  I  )
1739a04f6461SBarry Smith      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of
1740c8a0d604SMatthew G Knepley      ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given
1741c8a0d604SMatthew G Knepley      in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0,
1742c8a0d604SMatthew G Knepley      it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
1743a04f6461SBarry Smith      A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
1744c8a0d604SMatthew G Knepley      option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc. The
1745c9c6ffaaSJed Brown      factorization type is set using -pc_fieldsplit_schur_fact_type <diag, lower, upper, full>. The full is shown above,
17465668aaf4SBarry Smith      diag gives
1747c8a0d604SMatthew G Knepley $              ( inv(A00)     0   )
1748c8a0d604SMatthew G Knepley $              (   0      -ksp(S) )
17495668aaf4SBarry 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
1750c8a0d604SMatthew G Knepley $              (  A00   0 )
1751c8a0d604SMatthew G Knepley $              (  A10   S )
1752c8a0d604SMatthew G Knepley      where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of
1753c8a0d604SMatthew G Knepley $              ( A00 A01 )
1754c8a0d604SMatthew G Knepley $              (  0   S  )
1755c8a0d604SMatthew G Knepley      where again the inverses of A00 and S are applied using KSPs.
1756e69d4d44SBarry Smith 
1757edf189efSBarry Smith      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
1758edf189efSBarry Smith      is used automatically for a second block.
1759edf189efSBarry Smith 
1760ff218e97SBarry Smith      The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1.
1761ff218e97SBarry Smith      Generally it should be used with the AIJ format.
1762ff218e97SBarry Smith 
1763ff218e97SBarry Smith      The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see,
1764ff218e97SBarry Smith      for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT
1765ff218e97SBarry Smith      inside a smoother resulting in "Distributive Smoothers".
17660716a85fSBarry Smith 
1767a541d17aSBarry Smith    Concepts: physics based preconditioners, block preconditioners
17680971522cSBarry Smith 
17697e8cb189SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
1770e69d4d44SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
17710971522cSBarry Smith M*/
17720971522cSBarry Smith 
17730971522cSBarry Smith EXTERN_C_BEGIN
17740971522cSBarry Smith #undef __FUNCT__
17750971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
17767087cfbeSBarry Smith PetscErrorCode  PCCreate_FieldSplit(PC pc)
17770971522cSBarry Smith {
17780971522cSBarry Smith   PetscErrorCode ierr;
17790971522cSBarry Smith   PC_FieldSplit  *jac;
17800971522cSBarry Smith 
17810971522cSBarry Smith   PetscFunctionBegin;
178238f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
17832fa5cd67SKarl Rupp 
17840971522cSBarry Smith   jac->bs                 = -1;
17850971522cSBarry Smith   jac->nsplits            = 0;
17863e197d65SBarry Smith   jac->type               = PC_COMPOSITE_MULTIPLICATIVE;
1787e6cab6aaSJed Brown   jac->schurpre           = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
1788c9c6ffaaSJed Brown   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACT_FULL;
178951f519a2SBarry Smith 
17900971522cSBarry Smith   pc->data = (void*)jac;
17910971522cSBarry Smith 
17920971522cSBarry Smith   pc->ops->apply           = PCApply_FieldSplit;
1793421e10b8SBarry Smith   pc->ops->applytranspose  = PCApplyTranspose_FieldSplit;
17940971522cSBarry Smith   pc->ops->setup           = PCSetUp_FieldSplit;
1795574deadeSBarry Smith   pc->ops->reset           = PCReset_FieldSplit;
17960971522cSBarry Smith   pc->ops->destroy         = PCDestroy_FieldSplit;
17970971522cSBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_FieldSplit;
17980971522cSBarry Smith   pc->ops->view            = PCView_FieldSplit;
17990971522cSBarry Smith   pc->ops->applyrichardson = 0;
18000971522cSBarry Smith 
180169a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
180269a612a9SBarry Smith                                            PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
18030971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
18040971522cSBarry Smith                                            PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1805704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
1806704ba839SBarry Smith                                            PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
180779416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
180879416396SBarry Smith                                            PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
180951f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
181051f519a2SBarry Smith                                            PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
18110971522cSBarry Smith   PetscFunctionReturn(0);
18120971522cSBarry Smith }
18130971522cSBarry Smith EXTERN_C_END
18140971522cSBarry Smith 
18150971522cSBarry Smith 
1816a541d17aSBarry Smith 
1817