xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 5d4c12cdfef28496879fd13525d46a2a1b76f2dd)
1dba47a55SKris Buschelman 
2c6db04a5SJed Brown #include <private/pcimpl.h>     /*I "petscpc.h" I*/
33c48a1e8SJed Brown #include <petscdmcomposite.h>   /*I "petscdmcomposite.h" I*/
40971522cSBarry Smith 
5f7c28744SJed Brown const char *const PCFieldSplitSchurPreTypes[] = {"SELF","DIAG","USER","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0};
6f7c28744SJed Brown const char *const PCFieldSplitSchurFactorizationTypes[] = {"DIAG","LOWER","UPPER","FULL","PCFieldSplitSchurFactorizationType","PC_FIELDSPLIT_SCHUR_FACTORIZATION_",0};
7c5d2311dSJed Brown 
8c5d2311dSJed Brown typedef enum {
9c5d2311dSJed Brown   PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG,
10c5d2311dSJed Brown   PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER,
11c5d2311dSJed Brown   PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER,
12c5d2311dSJed Brown   PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL
13c5d2311dSJed Brown } PCFieldSplitSchurFactorizationType;
14084e4875SJed Brown 
150971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
160971522cSBarry Smith struct _PC_FieldSplitLink {
1769a612a9SBarry Smith   KSP               ksp;
180971522cSBarry Smith   Vec               x,y;
19db4c96c1SJed Brown   char              *splitname;
200971522cSBarry Smith   PetscInt          nfields;
21*5d4c12cdSJungho Lee   PetscInt          *fields,*fields_col;
221b9fc7fcSBarry Smith   VecScatter        sctx;
23*5d4c12cdSJungho Lee   IS                is,is_col;
2451f519a2SBarry Smith   PC_FieldSplitLink next,previous;
250971522cSBarry Smith };
260971522cSBarry Smith 
270971522cSBarry Smith typedef struct {
2868dd23aaSBarry Smith   PCCompositeType                    type;
29ace3abfcSBarry Smith   PetscBool                          defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */
30ace3abfcSBarry Smith   PetscBool                          splitdefined; /* Flag is set after the splits have been defined, to prevent more splits from being added */
31ace3abfcSBarry Smith   PetscBool                          realdiagonal; /* Flag to use the diagonal blocks of mat preconditioned by pmat, instead of just pmat */
3230ad9308SMatthew Knepley   PetscInt                           bs;           /* Block size for IS and Mat structures */
3330ad9308SMatthew Knepley   PetscInt                           nsplits;      /* Number of field divisions defined */
3479416396SBarry Smith   Vec                                *x,*y,w1,w2;
35519d70e2SJed Brown   Mat                                *mat;         /* The diagonal block for each split */
36519d70e2SJed Brown   Mat                                *pmat;        /* The preconditioning diagonal block for each split */
3730ad9308SMatthew Knepley   Mat                                *Afield;      /* The rows of the matrix associated with each split */
38ace3abfcSBarry Smith   PetscBool                          issetup;
3930ad9308SMatthew Knepley   /* Only used when Schur complement preconditioning is used */
4030ad9308SMatthew Knepley   Mat                                B;            /* The (0,1) block */
4130ad9308SMatthew Knepley   Mat                                C;            /* The (1,0) block */
42a04f6461SBarry Smith   Mat                                schur;        /* The Schur complement S = A11 - A10 A00^{-1} A01 */
43084e4875SJed Brown   Mat                                schur_user;   /* User-provided preconditioning matrix for the Schur complement */
44084e4875SJed Brown   PCFieldSplitSchurPreType           schurpre; /* Determines which preconditioning matrix is used for the Schur complement */
45c5d2311dSJed Brown   PCFieldSplitSchurFactorizationType schurfactorization;
4630ad9308SMatthew Knepley   KSP                                kspschur;     /* The solver for S */
4797bbdb24SBarry Smith   PC_FieldSplitLink                  head;
4863ec74ffSBarry Smith   PetscBool                          reset;         /* indicates PCReset() has been last called on this object, hack */
49c1570756SJed Brown   PetscBool                          suboptionsset; /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */
500971522cSBarry Smith } PC_FieldSplit;
510971522cSBarry Smith 
5216913363SBarry Smith /*
5316913363SBarry Smith     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
5416913363SBarry Smith    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
5516913363SBarry Smith    PC you could change this.
5616913363SBarry Smith */
57084e4875SJed Brown 
58e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it.  This way the
59084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */
60084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac)
61084e4875SJed Brown {
62084e4875SJed Brown   switch (jac->schurpre) {
63084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur;
64084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1];
65084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */
66084e4875SJed Brown     default:
67084e4875SJed Brown       return jac->schur_user ? jac->schur_user : jac->pmat[1];
68084e4875SJed Brown   }
69084e4875SJed Brown }
70084e4875SJed Brown 
71084e4875SJed Brown 
720971522cSBarry Smith #undef __FUNCT__
730971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
740971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
750971522cSBarry Smith {
760971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
770971522cSBarry Smith   PetscErrorCode    ierr;
78ace3abfcSBarry Smith   PetscBool         iascii;
790971522cSBarry Smith   PetscInt          i,j;
805a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
810971522cSBarry Smith 
820971522cSBarry Smith   PetscFunctionBegin;
832692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
840971522cSBarry Smith   if (iascii) {
8551f519a2SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
8669a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
870971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
880971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
891ab39975SBarry Smith       if (ilink->fields) {
900971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
9179416396SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
925a9f2f41SSatish Balay 	for (j=0; j<ilink->nfields; j++) {
9379416396SBarry Smith 	  if (j > 0) {
9479416396SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
9579416396SBarry Smith 	  }
965a9f2f41SSatish Balay 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
970971522cSBarry Smith 	}
980971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
9979416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1001ab39975SBarry Smith       } else {
1011ab39975SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1021ab39975SBarry Smith       }
1035a9f2f41SSatish Balay       ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
1045a9f2f41SSatish Balay       ilink = ilink->next;
1050971522cSBarry Smith     }
1060971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1070971522cSBarry Smith   } else {
10865e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1090971522cSBarry Smith   }
1100971522cSBarry Smith   PetscFunctionReturn(0);
1110971522cSBarry Smith }
1120971522cSBarry Smith 
1130971522cSBarry Smith #undef __FUNCT__
1143b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur"
1153b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
1163b224e63SBarry Smith {
1173b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1183b224e63SBarry Smith   PetscErrorCode    ierr;
119ace3abfcSBarry Smith   PetscBool         iascii;
1203b224e63SBarry Smith   PetscInt          i,j;
1213b224e63SBarry Smith   PC_FieldSplitLink ilink = jac->head;
1223b224e63SBarry Smith   KSP               ksp;
1233b224e63SBarry Smith 
1243b224e63SBarry Smith   PetscFunctionBegin;
1252692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1263b224e63SBarry Smith   if (iascii) {
127c5d2311dSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactorizationTypes[jac->schurfactorization]);CHKERRQ(ierr);
1283b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Split info:\n");CHKERRQ(ierr);
1293b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1303b224e63SBarry Smith     for (i=0; i<jac->nsplits; i++) {
1313b224e63SBarry Smith       if (ilink->fields) {
1323b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
1333b224e63SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
1343b224e63SBarry Smith 	for (j=0; j<ilink->nfields; j++) {
1353b224e63SBarry Smith 	  if (j > 0) {
1363b224e63SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
1373b224e63SBarry Smith 	  }
1383b224e63SBarry Smith 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
1393b224e63SBarry Smith 	}
1403b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1413b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1423b224e63SBarry Smith       } else {
1433b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1443b224e63SBarry Smith       }
1453b224e63SBarry Smith       ilink = ilink->next;
1463b224e63SBarry Smith     }
147435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block \n");CHKERRQ(ierr);
1483b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
14912cae6f2SJed Brown     if (jac->schur) {
15012cae6f2SJed Brown       ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
1513b224e63SBarry Smith       ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
15212cae6f2SJed Brown     } else {
15312cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
15412cae6f2SJed Brown     }
1553b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
156435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr);
1573b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15812cae6f2SJed Brown     if (jac->kspschur) {
1593b224e63SBarry Smith       ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
16012cae6f2SJed Brown     } else {
16112cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
16212cae6f2SJed Brown     }
1633b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1643b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1653b224e63SBarry Smith   } else {
16665e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1673b224e63SBarry Smith   }
1683b224e63SBarry Smith   PetscFunctionReturn(0);
1693b224e63SBarry Smith }
1703b224e63SBarry Smith 
1713b224e63SBarry Smith #undef __FUNCT__
1726c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private"
1736c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */
1746c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc)
1756c924f48SJed Brown {
1766c924f48SJed Brown   PetscErrorCode ierr;
1776c924f48SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
178*5d4c12cdSJungho Lee   PetscInt       i,nfields,*ifields,nfields_col,*ifields_col;
179*5d4c12cdSJungho Lee   PetscBool      flg,flg_col;
180*5d4c12cdSJungho Lee   char           optionname[128],splitname[8],optionname_col[128];
1816c924f48SJed Brown 
1826c924f48SJed Brown   PetscFunctionBegin;
1836c924f48SJed Brown   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr);
184*5d4c12cdSJungho Lee   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields_col);CHKERRQ(ierr);
1856c924f48SJed Brown   for (i=0,flg=PETSC_TRUE; ; i++) {
1866c924f48SJed Brown     ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
1876c924f48SJed Brown     ierr = PetscSNPrintf(optionname,sizeof optionname,"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr);
188*5d4c12cdSJungho Lee     ierr = PetscSNPrintf(optionname_col,sizeof optionname_col,"-pc_fieldsplit_%D_fields_col",i);CHKERRQ(ierr);
1896c924f48SJed Brown     nfields = jac->bs;
1906c924f48SJed Brown     ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr);
191*5d4c12cdSJungho Lee     ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname_col,ifields_col,&nfields_col,&flg_col);CHKERRQ(ierr);
1926c924f48SJed Brown     if (!flg) break;
193*5d4c12cdSJungho Lee     else if (flg && !flg_col) {
194*5d4c12cdSJungho Lee       if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
195*5d4c12cdSJungho Lee       ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields);CHKERRQ(ierr);
196*5d4c12cdSJungho Lee     }
197*5d4c12cdSJungho Lee     else {
198*5d4c12cdSJungho Lee       if (!nfields || !nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
199*5d4c12cdSJungho Lee       if (nfields != nfields_col) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Number of row and column fields must match");
200*5d4c12cdSJungho Lee       ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields,ifields_col);CHKERRQ(ierr);
201*5d4c12cdSJungho Lee     }
2026c924f48SJed Brown   }
2036c924f48SJed Brown   if (i > 0) {
2046c924f48SJed Brown     /* Makes command-line setting of splits take precedence over setting them in code.
2056c924f48SJed Brown        Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would
2066c924f48SJed Brown        create new splits, which would probably not be what the user wanted. */
2076c924f48SJed Brown     jac->splitdefined = PETSC_TRUE;
2086c924f48SJed Brown   }
2096c924f48SJed Brown   ierr = PetscFree(ifields);CHKERRQ(ierr);
210*5d4c12cdSJungho Lee   ierr = PetscFree(ifields_col);CHKERRQ(ierr);
2116c924f48SJed Brown   PetscFunctionReturn(0);
2126c924f48SJed Brown }
2136c924f48SJed Brown 
2146c924f48SJed Brown #undef __FUNCT__
21569a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
21669a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
2170971522cSBarry Smith {
2180971522cSBarry Smith   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
2190971522cSBarry Smith   PetscErrorCode    ierr;
2205a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
2216ce1633cSBarry Smith   PetscBool         flg = PETSC_FALSE,stokes = PETSC_FALSE;
2226c924f48SJed Brown   PetscInt          i;
2230971522cSBarry Smith 
2240971522cSBarry Smith   PetscFunctionBegin;
225d32f9abdSBarry Smith   if (!ilink) {
226d0af7cd3SBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
227d0af7cd3SBarry Smith     if (pc->dm && !stokes) {
2284d343eeaSMatthew G Knepley       PetscBool dmcomposite;
2294d343eeaSMatthew G Knepley       ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr);
2304d343eeaSMatthew G Knepley       if (dmcomposite) {
2314d343eeaSMatthew G Knepley         PetscInt nDM;
232*5d4c12cdSJungho Lee         IS       *fields;
2338b8307b2SJed Brown         ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
2348b8307b2SJed Brown         ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr);
235*5d4c12cdSJungho Lee         ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr);
236*5d4c12cdSJungho Lee         for (i=0; i<nDM; i++) {
237*5d4c12cdSJungho Lee           char splitname[8];
238*5d4c12cdSJungho Lee           ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
239*5d4c12cdSJungho Lee           ierr = PCFieldSplitSetIS(pc,splitname,fields[i]);CHKERRQ(ierr);
240*5d4c12cdSJungho Lee           ierr = ISDestroy(&fields[i]);CHKERRQ(ierr);
2412fa5ba8aSJed Brown         }
242*5d4c12cdSJungho Lee         ierr = PetscFree(fields);CHKERRQ(ierr);
2438b8307b2SJed Brown       }
24466ffff09SJed Brown     } else {
245521d7252SBarry Smith       if (jac->bs <= 0) {
246704ba839SBarry Smith         if (pc->pmat) {
247521d7252SBarry Smith           ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
248704ba839SBarry Smith         } else {
249704ba839SBarry Smith           jac->bs = 1;
250704ba839SBarry Smith         }
251521d7252SBarry Smith       }
252d32f9abdSBarry Smith 
253acfcf0e5SJed Brown       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
2546ce1633cSBarry Smith       if (stokes) {
2556ce1633cSBarry Smith         IS       zerodiags,rest;
2566ce1633cSBarry Smith         PetscInt nmin,nmax;
2576ce1633cSBarry Smith 
2586ce1633cSBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
2596ce1633cSBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
2606ce1633cSBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
2616ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
2626ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr);
263fcfd50ebSBarry Smith         ierr = ISDestroy(&zerodiags);CHKERRQ(ierr);
264fcfd50ebSBarry Smith         ierr = ISDestroy(&rest);CHKERRQ(ierr);
2656ce1633cSBarry Smith       } else {
266d32f9abdSBarry Smith         if (!flg) {
267d32f9abdSBarry Smith           /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
268d32f9abdSBarry Smith            then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
2696c924f48SJed Brown           ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
2706c924f48SJed Brown           if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
271d32f9abdSBarry Smith         }
2726c924f48SJed Brown         if (flg || !jac->splitdefined) {
273d32f9abdSBarry Smith           ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
274db4c96c1SJed Brown           for (i=0; i<jac->bs; i++) {
2756c924f48SJed Brown             char splitname[8];
2766c924f48SJed Brown             ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
277*5d4c12cdSJungho Lee             ierr = PCFieldSplitSetFields(pc,splitname,1,&i,&i);CHKERRQ(ierr);
27879416396SBarry Smith           }
279*5d4c12cdSJungho Lee           jac->defaultsplit = PETSC_TRUE;
280521d7252SBarry Smith         }
28166ffff09SJed Brown       }
2826ce1633cSBarry Smith     }
283edf189efSBarry Smith   } else if (jac->nsplits == 1) {
284edf189efSBarry Smith     if (ilink->is) {
285edf189efSBarry Smith       IS       is2;
286edf189efSBarry Smith       PetscInt nmin,nmax;
287edf189efSBarry Smith 
288edf189efSBarry Smith       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
289edf189efSBarry Smith       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
290db4c96c1SJed Brown       ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr);
291fcfd50ebSBarry Smith       ierr = ISDestroy(&is2);CHKERRQ(ierr);
292*5d4c12cdSJungho Lee     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
29363ec74ffSBarry Smith   } else if (jac->reset) {
29463ec74ffSBarry Smith     /* PCReset() has been called on this PC, ilink exists but all IS and Vec data structures in it must be rebuilt
295d0af7cd3SBarry Smith        This is basically the !ilink portion of code above copied from above and the allocation of the ilinks removed
296d0af7cd3SBarry Smith        since they already exist. This should be totally rewritten */
297d0af7cd3SBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
298d0af7cd3SBarry Smith     if (pc->dm && !stokes) {
299d0af7cd3SBarry Smith       PetscBool dmcomposite;
300d0af7cd3SBarry Smith       ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr);
301d0af7cd3SBarry Smith       if (dmcomposite) {
302d0af7cd3SBarry Smith         PetscInt nDM;
303d0af7cd3SBarry Smith         IS       *fields;
304d0af7cd3SBarry Smith         ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
305d0af7cd3SBarry Smith         ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr);
306d0af7cd3SBarry Smith         ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr);
307d0af7cd3SBarry Smith         for (i=0; i<nDM; i++) {
308d0af7cd3SBarry Smith           ilink->is = fields[i];
309d0af7cd3SBarry Smith           ilink     = ilink->next;
310edf189efSBarry Smith         }
311d0af7cd3SBarry Smith         ierr = PetscFree(fields);CHKERRQ(ierr);
312d0af7cd3SBarry Smith       }
313d0af7cd3SBarry Smith     } else {
314d0af7cd3SBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
315d0af7cd3SBarry Smith       if (stokes) {
316d0af7cd3SBarry Smith         IS       zerodiags,rest;
317d0af7cd3SBarry Smith         PetscInt nmin,nmax;
318d0af7cd3SBarry Smith 
319d0af7cd3SBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
320d0af7cd3SBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
321d0af7cd3SBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
32220b26d62SBarry Smith         ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
32320b26d62SBarry Smith         ierr = ISDestroy(&ilink->next->is);CHKERRQ(ierr);
324d0af7cd3SBarry Smith         ilink->is       = rest;
325d0af7cd3SBarry Smith         ilink->next->is = zerodiags;
32620b26d62SBarry Smith       } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used");
327d0af7cd3SBarry Smith     }
328d0af7cd3SBarry Smith   }
329d0af7cd3SBarry Smith 
330*5d4c12cdSJungho Lee   if (jac->nsplits < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields");
33169a612a9SBarry Smith   PetscFunctionReturn(0);
33269a612a9SBarry Smith }
33369a612a9SBarry Smith 
33469a612a9SBarry Smith #undef __FUNCT__
33569a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
33669a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
33769a612a9SBarry Smith {
33869a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
33969a612a9SBarry Smith   PetscErrorCode    ierr;
3405a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
341cf502942SBarry Smith   PetscInt          i,nsplit,ccsize;
34269a612a9SBarry Smith   MatStructure      flag = pc->flag;
343ace3abfcSBarry Smith   PetscBool         sorted;
34469a612a9SBarry Smith 
34569a612a9SBarry Smith   PetscFunctionBegin;
34669a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
34797bbdb24SBarry Smith   nsplit = jac->nsplits;
3485a9f2f41SSatish Balay   ilink  = jac->head;
34997bbdb24SBarry Smith 
35097bbdb24SBarry Smith   /* get the matrices for each split */
351704ba839SBarry Smith   if (!jac->issetup) {
3521b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
35397bbdb24SBarry Smith 
354704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
355704ba839SBarry Smith 
356*5d4c12cdSJungho Lee     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
35751f519a2SBarry Smith     bs     = jac->bs;
35897bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
359cf502942SBarry Smith     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
3601b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
3611b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
3621b9fc7fcSBarry Smith       if (jac->defaultsplit) {
363704ba839SBarry Smith         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
364704ba839SBarry Smith       } else if (!ilink->is) {
365ccb205f8SBarry Smith         if (ilink->nfields > 1) {
3665a9f2f41SSatish Balay           PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
3675a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
3681b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
3691b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
3701b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
37197bbdb24SBarry Smith             }
37297bbdb24SBarry Smith           }
373d67e408aSBarry Smith           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr);
374ccb205f8SBarry Smith         } else {
375704ba839SBarry Smith           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
376ccb205f8SBarry Smith         }
3773e197d65SBarry Smith       }
378704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
379e32f2f54SBarry Smith       if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
380704ba839SBarry Smith       ilink = ilink->next;
3811b9fc7fcSBarry Smith     }
3821b9fc7fcSBarry Smith   }
3831b9fc7fcSBarry Smith 
384704ba839SBarry Smith   ilink  = jac->head;
38597bbdb24SBarry Smith   if (!jac->pmat) {
386cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
387cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
3884aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
389704ba839SBarry Smith       ilink = ilink->next;
390cf502942SBarry Smith     }
39197bbdb24SBarry Smith   } else {
392cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
3934aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
394704ba839SBarry Smith       ilink = ilink->next;
395cf502942SBarry Smith     }
39697bbdb24SBarry Smith   }
397519d70e2SJed Brown   if (jac->realdiagonal) {
398519d70e2SJed Brown     ilink = jac->head;
399519d70e2SJed Brown     if (!jac->mat) {
400519d70e2SJed Brown       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr);
401519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
402519d70e2SJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
403519d70e2SJed Brown         ilink = ilink->next;
404519d70e2SJed Brown       }
405519d70e2SJed Brown     } else {
406519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
407966e74d7SJed Brown         if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);}
408519d70e2SJed Brown         ilink = ilink->next;
409519d70e2SJed Brown       }
410519d70e2SJed Brown     }
411519d70e2SJed Brown   } else {
412519d70e2SJed Brown     jac->mat = jac->pmat;
413519d70e2SJed Brown   }
41497bbdb24SBarry Smith 
4156c8605c2SJed Brown   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
41668dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
41768dd23aaSBarry Smith     ilink  = jac->head;
41868dd23aaSBarry Smith     if (!jac->Afield) {
41968dd23aaSBarry Smith       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
42068dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
4214aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
42268dd23aaSBarry Smith         ilink = ilink->next;
42368dd23aaSBarry Smith       }
42468dd23aaSBarry Smith     } else {
42568dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
4264aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
42768dd23aaSBarry Smith         ilink = ilink->next;
42868dd23aaSBarry Smith       }
42968dd23aaSBarry Smith     }
43068dd23aaSBarry Smith   }
43168dd23aaSBarry Smith 
4323b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
4333b224e63SBarry Smith     IS       ccis;
4344aa3045dSJed Brown     PetscInt rstart,rend;
435e7e72b3dSBarry Smith     if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
43668dd23aaSBarry Smith 
437e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
438e6cab6aaSJed Brown     ierr  = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
439e6cab6aaSJed Brown 
4403b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
4413b224e63SBarry Smith     if (jac->schur) {
4423b224e63SBarry Smith       ilink = jac->head;
4434aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4444aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
445fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4463b224e63SBarry Smith       ilink = ilink->next;
4474aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4484aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
449fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
450519d70e2SJed Brown       ierr  = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr);
451084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr);
4523b224e63SBarry Smith 
4533b224e63SBarry Smith      } else {
4541cee3971SBarry Smith       KSP ksp;
4556c924f48SJed Brown       char schurprefix[256];
4563b224e63SBarry Smith 
457a04f6461SBarry Smith       /* extract the A01 and A10 matrices */
4583b224e63SBarry Smith       ilink = jac->head;
4594aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4604aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
461fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4623b224e63SBarry Smith       ilink = ilink->next;
4634aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4644aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
465fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4667d50072eSJed Brown       /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */
467519d70e2SJed Brown       ierr  = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr);
468a04f6461SBarry Smith       /* set tabbing and options prefix of KSP inside the MatSchur */
4691cee3971SBarry Smith       ierr  = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
470e69d4d44SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr);
471a04f6461SBarry Smith       ierr  = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",jac->head->splitname);CHKERRQ(ierr);
472a04f6461SBarry Smith       ierr  = KSPSetOptionsPrefix(ksp,schurprefix);CHKERRQ(ierr);
47320b26d62SBarry Smith       /* Need to call this everytime because new matrix is being created */
47420b26d62SBarry Smith       ierr  = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
4753b224e63SBarry Smith 
4763b224e63SBarry Smith       ierr  = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr);
4779005cf84SBarry Smith       ierr  = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr);
4781cee3971SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
479084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
480084e4875SJed Brown       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
481084e4875SJed Brown         PC pc;
482cd5f4a64SJed Brown         ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr);
483084e4875SJed Brown         ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
484084e4875SJed Brown         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
485e69d4d44SBarry Smith       }
4866c924f48SJed Brown       ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
4876c924f48SJed Brown       ierr  = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr);
4883b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
48920b26d62SBarry Smith       /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */
49020b26d62SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
4913b224e63SBarry Smith 
4923b224e63SBarry Smith       ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr);
4933b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr);
4943b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr);
4953b224e63SBarry Smith       ilink = jac->head;
4963b224e63SBarry Smith       ilink->x = jac->x[0]; ilink->y = jac->y[0];
4973b224e63SBarry Smith       ilink = ilink->next;
4983b224e63SBarry Smith       ilink->x = jac->x[1]; ilink->y = jac->y[1];
4993b224e63SBarry Smith     }
5003b224e63SBarry Smith   } else {
50197bbdb24SBarry Smith     /* set up the individual PCs */
50297bbdb24SBarry Smith     i    = 0;
5035a9f2f41SSatish Balay     ilink = jac->head;
5045a9f2f41SSatish Balay     while (ilink) {
505519d70e2SJed Brown       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr);
5063b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
507c1570756SJed Brown       if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);}
5085a9f2f41SSatish Balay       ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
50997bbdb24SBarry Smith       i++;
5105a9f2f41SSatish Balay       ilink = ilink->next;
5110971522cSBarry Smith     }
51297bbdb24SBarry Smith 
51397bbdb24SBarry Smith     /* create work vectors for each split */
5141b9fc7fcSBarry Smith     if (!jac->x) {
51597bbdb24SBarry Smith       ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
5165a9f2f41SSatish Balay       ilink = jac->head;
51797bbdb24SBarry Smith       for (i=0; i<nsplit; i++) {
518906ed7ccSBarry Smith         Vec *vl,*vr;
519906ed7ccSBarry Smith 
520906ed7ccSBarry Smith         ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
521906ed7ccSBarry Smith         ilink->x  = *vr;
522906ed7ccSBarry Smith         ilink->y  = *vl;
523906ed7ccSBarry Smith         ierr      = PetscFree(vr);CHKERRQ(ierr);
524906ed7ccSBarry Smith         ierr      = PetscFree(vl);CHKERRQ(ierr);
5255a9f2f41SSatish Balay         jac->x[i] = ilink->x;
5265a9f2f41SSatish Balay         jac->y[i] = ilink->y;
5275a9f2f41SSatish Balay         ilink     = ilink->next;
52897bbdb24SBarry Smith       }
5293b224e63SBarry Smith     }
5303b224e63SBarry Smith   }
5313b224e63SBarry Smith 
5323b224e63SBarry Smith 
533d0f46423SBarry Smith   if (!jac->head->sctx) {
5343b224e63SBarry Smith     Vec xtmp;
5353b224e63SBarry Smith 
53679416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
5371b9fc7fcSBarry Smith 
5385a9f2f41SSatish Balay     ilink = jac->head;
5391b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
5401b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
541704ba839SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
5425a9f2f41SSatish Balay       ilink = ilink->next;
5431b9fc7fcSBarry Smith     }
5446bf464f9SBarry Smith     ierr = VecDestroy(&xtmp);CHKERRQ(ierr);
5451b9fc7fcSBarry Smith   }
546c1570756SJed Brown   jac->suboptionsset = PETSC_TRUE;
5470971522cSBarry Smith   PetscFunctionReturn(0);
5480971522cSBarry Smith }
5490971522cSBarry Smith 
5505a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
551ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
552ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
5535a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
554ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
555ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
55679416396SBarry Smith 
5570971522cSBarry Smith #undef __FUNCT__
5583b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur"
5593b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
5603b224e63SBarry Smith {
5613b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
5623b224e63SBarry Smith   PetscErrorCode    ierr;
5633b224e63SBarry Smith   KSP               ksp;
5643b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
5653b224e63SBarry Smith 
5663b224e63SBarry Smith   PetscFunctionBegin;
5673b224e63SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
5683b224e63SBarry Smith 
569c5d2311dSJed Brown   switch (jac->schurfactorization) {
570c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG:
571a04f6461SBarry Smith       /* [A00 0; 0 -S], positive definite, suitable for MINRES */
572c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
573c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
574c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
575c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
576c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
577c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
578c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
579c5d2311dSJed Brown       ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr);
580c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
581c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
582c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
583c5d2311dSJed Brown       break;
584c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER:
585a04f6461SBarry Smith       /* [A00 0; A10 S], suitable for left preconditioning */
586c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
587c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
588c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
589c5d2311dSJed Brown       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
590c5d2311dSJed Brown       ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr);
591c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
592c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
593c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
594c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
595c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
596c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
597c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
598c5d2311dSJed Brown       break;
599c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER:
600a04f6461SBarry Smith       /* [A00 A01; 0 S], suitable for right preconditioning */
601c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
602c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
603c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
604c5d2311dSJed Brown       ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
605c5d2311dSJed Brown       ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr);
606c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
607c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
608c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
609c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
610c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
611c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
612c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
613c5d2311dSJed Brown       break;
614c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL:
615a04f6461SBarry 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 */
6163b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6173b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6183b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
6193b224e63SBarry Smith       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
6203b224e63SBarry Smith       ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
6213b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6223b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6233b224e63SBarry Smith 
6243b224e63SBarry Smith       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
6253b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6263b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6273b224e63SBarry Smith 
6283b224e63SBarry Smith       ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
6293b224e63SBarry Smith       ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
6303b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
6313b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6323b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
633c5d2311dSJed Brown   }
6343b224e63SBarry Smith   PetscFunctionReturn(0);
6353b224e63SBarry Smith }
6363b224e63SBarry Smith 
6373b224e63SBarry Smith #undef __FUNCT__
6380971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
6390971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
6400971522cSBarry Smith {
6410971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
6420971522cSBarry Smith   PetscErrorCode    ierr;
6435a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
644af93645bSJed Brown   PetscInt          cnt;
6450971522cSBarry Smith 
6460971522cSBarry Smith   PetscFunctionBegin;
64751f519a2SBarry Smith   CHKMEMQ;
64851f519a2SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
64951f519a2SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
65051f519a2SBarry Smith 
65179416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
6521b9fc7fcSBarry Smith     if (jac->defaultsplit) {
6530971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
6545a9f2f41SSatish Balay       while (ilink) {
6555a9f2f41SSatish Balay         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
6565a9f2f41SSatish Balay         ilink = ilink->next;
6570971522cSBarry Smith       }
6580971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
6591b9fc7fcSBarry Smith     } else {
660efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
6615a9f2f41SSatish Balay       while (ilink) {
6625a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
6635a9f2f41SSatish Balay         ilink = ilink->next;
6641b9fc7fcSBarry Smith       }
6651b9fc7fcSBarry Smith     }
66616913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
66779416396SBarry Smith     if (!jac->w1) {
66879416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
66979416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
67079416396SBarry Smith     }
671efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
6725a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
6733e197d65SBarry Smith     cnt = 1;
6745a9f2f41SSatish Balay     while (ilink->next) {
6755a9f2f41SSatish Balay       ilink = ilink->next;
6763e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
6773e197d65SBarry Smith       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
6783e197d65SBarry Smith       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
6793e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6803e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6813e197d65SBarry Smith       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
6823e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6833e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6843e197d65SBarry Smith     }
68551f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
68611755939SBarry Smith       cnt -= 2;
68751f519a2SBarry Smith       while (ilink->previous) {
68851f519a2SBarry Smith         ilink = ilink->previous;
68911755939SBarry Smith         /* compute the residual only over the part of the vector needed */
69011755939SBarry Smith         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
69111755939SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
69211755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
69311755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
69411755939SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
69511755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
69611755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
69751f519a2SBarry Smith       }
69811755939SBarry Smith     }
69965e19b50SBarry Smith   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
70051f519a2SBarry Smith   CHKMEMQ;
7010971522cSBarry Smith   PetscFunctionReturn(0);
7020971522cSBarry Smith }
7030971522cSBarry Smith 
704421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
705ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
706ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
707421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
708ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
709ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
710421e10b8SBarry Smith 
711421e10b8SBarry Smith #undef __FUNCT__
7128c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit"
713421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
714421e10b8SBarry Smith {
715421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
716421e10b8SBarry Smith   PetscErrorCode    ierr;
717421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
718421e10b8SBarry Smith 
719421e10b8SBarry Smith   PetscFunctionBegin;
720421e10b8SBarry Smith   CHKMEMQ;
721421e10b8SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
722421e10b8SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
723421e10b8SBarry Smith 
724421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
725421e10b8SBarry Smith     if (jac->defaultsplit) {
726421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
727421e10b8SBarry Smith       while (ilink) {
728421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
729421e10b8SBarry Smith 	ilink = ilink->next;
730421e10b8SBarry Smith       }
731421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
732421e10b8SBarry Smith     } else {
733421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
734421e10b8SBarry Smith       while (ilink) {
735421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
736421e10b8SBarry Smith 	ilink = ilink->next;
737421e10b8SBarry Smith       }
738421e10b8SBarry Smith     }
739421e10b8SBarry Smith   } else {
740421e10b8SBarry Smith     if (!jac->w1) {
741421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
742421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
743421e10b8SBarry Smith     }
744421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
745421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
746421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
747421e10b8SBarry Smith       while (ilink->next) {
748421e10b8SBarry Smith         ilink = ilink->next;
7499989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
750421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
751421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
752421e10b8SBarry Smith       }
753421e10b8SBarry Smith       while (ilink->previous) {
754421e10b8SBarry Smith         ilink = ilink->previous;
7559989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
756421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
757421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
758421e10b8SBarry Smith       }
759421e10b8SBarry Smith     } else {
760421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
761421e10b8SBarry Smith 	ilink = ilink->next;
762421e10b8SBarry Smith       }
763421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
764421e10b8SBarry Smith       while (ilink->previous) {
765421e10b8SBarry Smith 	ilink = ilink->previous;
7669989ab13SBarry Smith 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
767421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
768421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
769421e10b8SBarry Smith       }
770421e10b8SBarry Smith     }
771421e10b8SBarry Smith   }
772421e10b8SBarry Smith   CHKMEMQ;
773421e10b8SBarry Smith   PetscFunctionReturn(0);
774421e10b8SBarry Smith }
775421e10b8SBarry Smith 
7760971522cSBarry Smith #undef __FUNCT__
777574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit"
778574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc)
7790971522cSBarry Smith {
7800971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
7810971522cSBarry Smith   PetscErrorCode    ierr;
7825a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
7830971522cSBarry Smith 
7840971522cSBarry Smith   PetscFunctionBegin;
7855a9f2f41SSatish Balay   while (ilink) {
786574deadeSBarry Smith     ierr = KSPReset(ilink->ksp);CHKERRQ(ierr);
787fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->x);CHKERRQ(ierr);
788fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->y);CHKERRQ(ierr);
789fcfd50ebSBarry Smith     ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr);
790fcfd50ebSBarry Smith     ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
7915a9f2f41SSatish Balay     next = ilink->next;
7925a9f2f41SSatish Balay     ilink = next;
7930971522cSBarry Smith   }
79405b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
795574deadeSBarry Smith   if (jac->mat && jac->mat != jac->pmat) {
796574deadeSBarry Smith     ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);
797574deadeSBarry Smith   } else if (jac->mat) {
798574deadeSBarry Smith     jac->mat = PETSC_NULL;
799574deadeSBarry Smith   }
80097bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
80168dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
8026bf464f9SBarry Smith   ierr = VecDestroy(&jac->w1);CHKERRQ(ierr);
8036bf464f9SBarry Smith   ierr = VecDestroy(&jac->w2);CHKERRQ(ierr);
8046bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur);CHKERRQ(ierr);
8056bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
8066bf464f9SBarry Smith   ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr);
8076bf464f9SBarry Smith   ierr = MatDestroy(&jac->B);CHKERRQ(ierr);
8086bf464f9SBarry Smith   ierr = MatDestroy(&jac->C);CHKERRQ(ierr);
80963ec74ffSBarry Smith   jac->reset = PETSC_TRUE;
810574deadeSBarry Smith   PetscFunctionReturn(0);
811574deadeSBarry Smith }
812574deadeSBarry Smith 
813574deadeSBarry Smith #undef __FUNCT__
814574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
815574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
816574deadeSBarry Smith {
817574deadeSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
818574deadeSBarry Smith   PetscErrorCode    ierr;
819574deadeSBarry Smith   PC_FieldSplitLink ilink = jac->head,next;
820574deadeSBarry Smith 
821574deadeSBarry Smith   PetscFunctionBegin;
822574deadeSBarry Smith   ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr);
823574deadeSBarry Smith   while (ilink) {
8246bf464f9SBarry Smith     ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr);
825574deadeSBarry Smith     next = ilink->next;
826574deadeSBarry Smith     ierr = PetscFree(ilink->splitname);CHKERRQ(ierr);
827574deadeSBarry Smith     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
828*5d4c12cdSJungho Lee     ierr = PetscFree(ilink->fields_col);CHKERRQ(ierr);
829574deadeSBarry Smith     ierr = PetscFree(ilink);CHKERRQ(ierr);
830574deadeSBarry Smith     ilink = next;
831574deadeSBarry Smith   }
832574deadeSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
833c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
8340971522cSBarry Smith   PetscFunctionReturn(0);
8350971522cSBarry Smith }
8360971522cSBarry Smith 
8370971522cSBarry Smith #undef __FUNCT__
8380971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
8390971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
8400971522cSBarry Smith {
8411b9fc7fcSBarry Smith   PetscErrorCode  ierr;
8426c924f48SJed Brown   PetscInt        bs;
843bc59fbc5SBarry Smith   PetscBool       flg,stokes = PETSC_FALSE;
8449dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
8453b224e63SBarry Smith   PCCompositeType ctype;
8461b9fc7fcSBarry Smith 
8470971522cSBarry Smith   PetscFunctionBegin;
8481b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
849acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr);
85051f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
85151f519a2SBarry Smith   if (flg) {
85251f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
85351f519a2SBarry Smith   }
854704ba839SBarry Smith 
855435f959eSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
856c0adfefeSBarry Smith   if (stokes) {
857c0adfefeSBarry Smith     ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
858c0adfefeSBarry Smith     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
859c0adfefeSBarry Smith   }
860c0adfefeSBarry Smith 
8613b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
8623b224e63SBarry Smith   if (flg) {
8633b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
8643b224e63SBarry Smith   }
865d32f9abdSBarry Smith 
866c30613efSMatthew Knepley   /* Only setup fields once */
867c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
868d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
869d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
8706c924f48SJed Brown     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
8716c924f48SJed Brown     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
872d32f9abdSBarry Smith   }
873c5d2311dSJed Brown   if (jac->type == PC_COMPOSITE_SCHUR) {
874c5d2311dSJed Brown     ierr = PetscOptionsEnum("-pc_fieldsplit_schur_factorization_type","Factorization to use","None",PCFieldSplitSchurFactorizationTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,PETSC_NULL);CHKERRQ(ierr);
875084e4875SJed Brown     ierr = PetscOptionsEnum("-pc_fieldsplit_schur_precondition","How to build preconditioner for Schur complement","PCFieldSplitSchurPrecondition",PCFieldSplitSchurPreTypes,(PetscEnum)jac->schurpre,(PetscEnum*)&jac->schurpre,PETSC_NULL);CHKERRQ(ierr);
876c5d2311dSJed Brown   }
8771b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
8780971522cSBarry Smith   PetscFunctionReturn(0);
8790971522cSBarry Smith }
8800971522cSBarry Smith 
8810971522cSBarry Smith /*------------------------------------------------------------------------------------*/
8820971522cSBarry Smith 
8830971522cSBarry Smith EXTERN_C_BEGIN
8840971522cSBarry Smith #undef __FUNCT__
8850971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
886*5d4c12cdSJungho Lee PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
8870971522cSBarry Smith {
88897bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
8890971522cSBarry Smith   PetscErrorCode    ierr;
8905a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
89169a612a9SBarry Smith   char              prefix[128];
892*5d4c12cdSJungho Lee   PetscInt          i;
8930971522cSBarry Smith 
8940971522cSBarry Smith   PetscFunctionBegin;
8956c924f48SJed Brown   if (jac->splitdefined) {
8966c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
8976c924f48SJed Brown     PetscFunctionReturn(0);
8986c924f48SJed Brown   }
89951f519a2SBarry Smith   for (i=0; i<n; i++) {
900e32f2f54SBarry 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);
901e32f2f54SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
90251f519a2SBarry Smith   }
903704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
904a04f6461SBarry Smith   if (splitname) {
905db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
906a04f6461SBarry Smith   } else {
907a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
908a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
909a04f6461SBarry Smith   }
910704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
9115a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
912*5d4c12cdSJungho Lee   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields_col);CHKERRQ(ierr);
913*5d4c12cdSJungho Lee   ierr = PetscMemcpy(ilink->fields_col,fields_col,n*sizeof(PetscInt));CHKERRQ(ierr);
9145a9f2f41SSatish Balay   ilink->nfields = n;
9155a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
9167adad957SLisandro Dalcin   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
9171cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
9185a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
9199005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
92069a612a9SBarry Smith 
921a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
9225a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
9230971522cSBarry Smith 
9240971522cSBarry Smith   if (!next) {
9255a9f2f41SSatish Balay     jac->head       = ilink;
92651f519a2SBarry Smith     ilink->previous = PETSC_NULL;
9270971522cSBarry Smith   } else {
9280971522cSBarry Smith     while (next->next) {
9290971522cSBarry Smith       next = next->next;
9300971522cSBarry Smith     }
9315a9f2f41SSatish Balay     next->next      = ilink;
93251f519a2SBarry Smith     ilink->previous = next;
9330971522cSBarry Smith   }
9340971522cSBarry Smith   jac->nsplits++;
9350971522cSBarry Smith   PetscFunctionReturn(0);
9360971522cSBarry Smith }
9370971522cSBarry Smith EXTERN_C_END
9380971522cSBarry Smith 
939e69d4d44SBarry Smith EXTERN_C_BEGIN
940e69d4d44SBarry Smith #undef __FUNCT__
941e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
9427087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
943e69d4d44SBarry Smith {
944e69d4d44SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
945e69d4d44SBarry Smith   PetscErrorCode ierr;
946e69d4d44SBarry Smith 
947e69d4d44SBarry Smith   PetscFunctionBegin;
948519d70e2SJed Brown   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
949e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
950e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
95113e0d083SBarry Smith   if (n) *n = jac->nsplits;
952e69d4d44SBarry Smith   PetscFunctionReturn(0);
953e69d4d44SBarry Smith }
954e69d4d44SBarry Smith EXTERN_C_END
9550971522cSBarry Smith 
9560971522cSBarry Smith EXTERN_C_BEGIN
9570971522cSBarry Smith #undef __FUNCT__
95869a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
9597087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
9600971522cSBarry Smith {
9610971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
9620971522cSBarry Smith   PetscErrorCode    ierr;
9630971522cSBarry Smith   PetscInt          cnt = 0;
9645a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
9650971522cSBarry Smith 
9660971522cSBarry Smith   PetscFunctionBegin;
9675d480477SMatthew G Knepley   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
9685a9f2f41SSatish Balay   while (ilink) {
9695a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
9705a9f2f41SSatish Balay     ilink = ilink->next;
9710971522cSBarry Smith   }
9725d480477SMatthew 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);
97313e0d083SBarry Smith   if (n) *n = jac->nsplits;
9740971522cSBarry Smith   PetscFunctionReturn(0);
9750971522cSBarry Smith }
9760971522cSBarry Smith EXTERN_C_END
9770971522cSBarry Smith 
978704ba839SBarry Smith EXTERN_C_BEGIN
979704ba839SBarry Smith #undef __FUNCT__
980704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
9817087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
982704ba839SBarry Smith {
983704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
984704ba839SBarry Smith   PetscErrorCode    ierr;
985704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
986704ba839SBarry Smith   char              prefix[128];
987704ba839SBarry Smith 
988704ba839SBarry Smith   PetscFunctionBegin;
9896c924f48SJed Brown   if (jac->splitdefined) {
9906c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
9916c924f48SJed Brown     PetscFunctionReturn(0);
9926c924f48SJed Brown   }
99316913363SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
994a04f6461SBarry Smith   if (splitname) {
995db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
996a04f6461SBarry Smith   } else {
997a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
998*5d4c12cdSJungho Lee     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
999a04f6461SBarry Smith   }
10001ab39975SBarry Smith   ilink->is      = is;
10011ab39975SBarry Smith   ierr           = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1002704ba839SBarry Smith   ilink->next    = PETSC_NULL;
1003704ba839SBarry Smith   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
10041cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1005704ba839SBarry Smith   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
10069005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
1007704ba839SBarry Smith 
1008a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
1009704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
1010704ba839SBarry Smith 
1011704ba839SBarry Smith   if (!next) {
1012704ba839SBarry Smith     jac->head       = ilink;
1013704ba839SBarry Smith     ilink->previous = PETSC_NULL;
1014704ba839SBarry Smith   } else {
1015704ba839SBarry Smith     while (next->next) {
1016704ba839SBarry Smith       next = next->next;
1017704ba839SBarry Smith     }
1018704ba839SBarry Smith     next->next      = ilink;
1019704ba839SBarry Smith     ilink->previous = next;
1020704ba839SBarry Smith   }
1021704ba839SBarry Smith   jac->nsplits++;
1022704ba839SBarry Smith 
1023704ba839SBarry Smith   PetscFunctionReturn(0);
1024704ba839SBarry Smith }
1025704ba839SBarry Smith EXTERN_C_END
1026704ba839SBarry Smith 
10270971522cSBarry Smith #undef __FUNCT__
10280971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
10290971522cSBarry Smith /*@
10300971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
10310971522cSBarry Smith 
1032ad4df100SBarry Smith     Logically Collective on PC
10330971522cSBarry Smith 
10340971522cSBarry Smith     Input Parameters:
10350971522cSBarry Smith +   pc  - the preconditioner context
1036a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
10370971522cSBarry Smith .   n - the number of fields in this split
1038db4c96c1SJed Brown -   fields - the fields in this split
10390971522cSBarry Smith 
10400971522cSBarry Smith     Level: intermediate
10410971522cSBarry Smith 
1042d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
1043d32f9abdSBarry Smith 
1044d32f9abdSBarry Smith      The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block
1045d32f9abdSBarry 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
1046d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1047d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
1048d32f9abdSBarry Smith 
1049db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
1050db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1051db4c96c1SJed Brown 
1052*5d4c12cdSJungho Lee      Developer Note: This routine does not actually create the IS representing the split, that is delayed
1053*5d4c12cdSJungho Lee      until PCSetUp_FieldSplit(), because information about the vector/matrix layouts may not be
1054*5d4c12cdSJungho Lee      available when this routine is called.
1055*5d4c12cdSJungho Lee 
1056d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
10570971522cSBarry Smith 
10580971522cSBarry Smith @*/
1059*5d4c12cdSJungho Lee PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields,const PetscInt *fields_col)
10600971522cSBarry Smith {
10614ac538c5SBarry Smith   PetscErrorCode ierr;
10620971522cSBarry Smith 
10630971522cSBarry Smith   PetscFunctionBegin;
10640700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1065db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1066db4c96c1SJed Brown   if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1067db4c96c1SJed Brown   PetscValidIntPointer(fields,3);
1068*5d4c12cdSJungho Lee   ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *,const PetscInt *),(pc,splitname,n,fields,fields_col));CHKERRQ(ierr);
10690971522cSBarry Smith   PetscFunctionReturn(0);
10700971522cSBarry Smith }
10710971522cSBarry Smith 
10720971522cSBarry Smith #undef __FUNCT__
1073704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
1074704ba839SBarry Smith /*@
1075704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
1076704ba839SBarry Smith 
1077ad4df100SBarry Smith     Logically Collective on PC
1078704ba839SBarry Smith 
1079704ba839SBarry Smith     Input Parameters:
1080704ba839SBarry Smith +   pc  - the preconditioner context
1081a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
1082db4c96c1SJed Brown -   is - the index set that defines the vector elements in this field
1083704ba839SBarry Smith 
1084d32f9abdSBarry Smith 
1085a6ffb8dbSJed Brown     Notes:
1086a6ffb8dbSJed Brown     Use PCFieldSplitSetFields(), for fields defined by strided types.
1087a6ffb8dbSJed Brown 
1088db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
1089db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1090d32f9abdSBarry Smith 
1091704ba839SBarry Smith     Level: intermediate
1092704ba839SBarry Smith 
1093704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
1094704ba839SBarry Smith 
1095704ba839SBarry Smith @*/
10967087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1097704ba839SBarry Smith {
10984ac538c5SBarry Smith   PetscErrorCode ierr;
1099704ba839SBarry Smith 
1100704ba839SBarry Smith   PetscFunctionBegin;
11010700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1102*5d4c12cdSJungho Lee   PetscValidCharPointer(splitname,2);
1103db4c96c1SJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,3);
11044ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr);
1105704ba839SBarry Smith   PetscFunctionReturn(0);
1106704ba839SBarry Smith }
1107704ba839SBarry Smith 
1108704ba839SBarry Smith #undef __FUNCT__
110957a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS"
111057a9adfeSMatthew G Knepley /*@
111157a9adfeSMatthew G Knepley     PCFieldSplitGetIS - Retrieves the elements for a field as an IS
111257a9adfeSMatthew G Knepley 
111357a9adfeSMatthew G Knepley     Logically Collective on PC
111457a9adfeSMatthew G Knepley 
111557a9adfeSMatthew G Knepley     Input Parameters:
111657a9adfeSMatthew G Knepley +   pc  - the preconditioner context
111757a9adfeSMatthew G Knepley -   splitname - name of this split
111857a9adfeSMatthew G Knepley 
111957a9adfeSMatthew G Knepley     Output Parameter:
112057a9adfeSMatthew G Knepley -   is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found
112157a9adfeSMatthew G Knepley 
112257a9adfeSMatthew G Knepley     Level: intermediate
112357a9adfeSMatthew G Knepley 
112457a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS()
112557a9adfeSMatthew G Knepley 
112657a9adfeSMatthew G Knepley @*/
112757a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
112857a9adfeSMatthew G Knepley {
112957a9adfeSMatthew G Knepley   PetscErrorCode ierr;
113057a9adfeSMatthew G Knepley 
113157a9adfeSMatthew G Knepley   PetscFunctionBegin;
113257a9adfeSMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
113357a9adfeSMatthew G Knepley   PetscValidCharPointer(splitname,2);
113457a9adfeSMatthew G Knepley   PetscValidPointer(is,3);
113557a9adfeSMatthew G Knepley   {
113657a9adfeSMatthew G Knepley     PC_FieldSplit    *jac   = (PC_FieldSplit *) pc->data;
113757a9adfeSMatthew G Knepley     PC_FieldSplitLink ilink = jac->head;
113857a9adfeSMatthew G Knepley     PetscBool         found;
113957a9adfeSMatthew G Knepley 
114057a9adfeSMatthew G Knepley     *is = PETSC_NULL;
114157a9adfeSMatthew G Knepley     while(ilink) {
114257a9adfeSMatthew G Knepley       ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr);
114357a9adfeSMatthew G Knepley       if (found) {
114457a9adfeSMatthew G Knepley         *is = ilink->is;
114557a9adfeSMatthew G Knepley         break;
114657a9adfeSMatthew G Knepley       }
114757a9adfeSMatthew G Knepley       ilink = ilink->next;
114857a9adfeSMatthew G Knepley     }
114957a9adfeSMatthew G Knepley   }
115057a9adfeSMatthew G Knepley   PetscFunctionReturn(0);
115157a9adfeSMatthew G Knepley }
115257a9adfeSMatthew G Knepley 
115357a9adfeSMatthew G Knepley #undef __FUNCT__
115451f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
115551f519a2SBarry Smith /*@
115651f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
115751f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
115851f519a2SBarry Smith 
1159ad4df100SBarry Smith     Logically Collective on PC
116051f519a2SBarry Smith 
116151f519a2SBarry Smith     Input Parameters:
116251f519a2SBarry Smith +   pc  - the preconditioner context
116351f519a2SBarry Smith -   bs - the block size
116451f519a2SBarry Smith 
116551f519a2SBarry Smith     Level: intermediate
116651f519a2SBarry Smith 
116751f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
116851f519a2SBarry Smith 
116951f519a2SBarry Smith @*/
11707087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
117151f519a2SBarry Smith {
11724ac538c5SBarry Smith   PetscErrorCode ierr;
117351f519a2SBarry Smith 
117451f519a2SBarry Smith   PetscFunctionBegin;
11750700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1176c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,bs,2);
11774ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr);
117851f519a2SBarry Smith   PetscFunctionReturn(0);
117951f519a2SBarry Smith }
118051f519a2SBarry Smith 
118151f519a2SBarry Smith #undef __FUNCT__
118269a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
11830971522cSBarry Smith /*@C
118469a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
11850971522cSBarry Smith 
118669a612a9SBarry Smith    Collective on KSP
11870971522cSBarry Smith 
11880971522cSBarry Smith    Input Parameter:
11890971522cSBarry Smith .  pc - the preconditioner context
11900971522cSBarry Smith 
11910971522cSBarry Smith    Output Parameters:
119213e0d083SBarry Smith +  n - the number of splits
119369a612a9SBarry Smith -  pc - the array of KSP contexts
11940971522cSBarry Smith 
11950971522cSBarry Smith    Note:
1196d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
1197d32f9abdSBarry Smith    (not the KSP just the array that contains them).
11980971522cSBarry Smith 
119969a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
12000971522cSBarry Smith 
12010971522cSBarry Smith    Level: advanced
12020971522cSBarry Smith 
12030971522cSBarry Smith .seealso: PCFIELDSPLIT
12040971522cSBarry Smith @*/
12057087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
12060971522cSBarry Smith {
12074ac538c5SBarry Smith   PetscErrorCode ierr;
12080971522cSBarry Smith 
12090971522cSBarry Smith   PetscFunctionBegin;
12100700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
121113e0d083SBarry Smith   if (n) PetscValidIntPointer(n,2);
12124ac538c5SBarry Smith   ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
12130971522cSBarry Smith   PetscFunctionReturn(0);
12140971522cSBarry Smith }
12150971522cSBarry Smith 
1216e69d4d44SBarry Smith #undef __FUNCT__
1217e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition"
1218e69d4d44SBarry Smith /*@
1219e69d4d44SBarry Smith     PCFieldSplitSchurPrecondition -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
1220a04f6461SBarry Smith       A11 matrix. Otherwise no preconditioner is used.
1221e69d4d44SBarry Smith 
1222e69d4d44SBarry Smith     Collective on PC
1223e69d4d44SBarry Smith 
1224e69d4d44SBarry Smith     Input Parameters:
1225e69d4d44SBarry Smith +   pc  - the preconditioner context
1226fd1303e9SJungho Lee .   ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default
1227084e4875SJed Brown -   userpre - matrix to use for preconditioning, or PETSC_NULL
1228084e4875SJed Brown 
1229e69d4d44SBarry Smith     Options Database:
1230084e4875SJed Brown .     -pc_fieldsplit_schur_precondition <self,user,diag> default is diag
1231e69d4d44SBarry Smith 
1232fd1303e9SJungho Lee     Notes:
1233fd1303e9SJungho Lee $    If ptype is
1234fd1303e9SJungho Lee $        user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument
1235fd1303e9SJungho Lee $             to this function).
1236fd1303e9SJungho Lee $        diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original
1237fd1303e9SJungho Lee $             matrix associated with the Schur complement (i.e. A11)
1238fd1303e9SJungho Lee $        self the preconditioner for the Schur complement is generated from the Schur complement matrix itself:
1239fd1303e9SJungho Lee $             The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC
1240fd1303e9SJungho Lee $             preconditioner
1241fd1303e9SJungho Lee 
1242fd1303e9SJungho Lee      When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense
1243fd1303e9SJungho Lee     with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
1244fd1303e9SJungho Lee     -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement.
1245fd1303e9SJungho Lee 
1246fd1303e9SJungho Lee     Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in
1247fd1303e9SJungho Lee      the name since it sets a proceedure to use.
1248fd1303e9SJungho Lee 
1249e69d4d44SBarry Smith     Level: intermediate
1250e69d4d44SBarry Smith 
1251fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC
1252e69d4d44SBarry Smith 
1253e69d4d44SBarry Smith @*/
12547087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1255e69d4d44SBarry Smith {
12564ac538c5SBarry Smith   PetscErrorCode ierr;
1257e69d4d44SBarry Smith 
1258e69d4d44SBarry Smith   PetscFunctionBegin;
12590700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12604ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr);
1261e69d4d44SBarry Smith   PetscFunctionReturn(0);
1262e69d4d44SBarry Smith }
1263e69d4d44SBarry Smith 
1264e69d4d44SBarry Smith EXTERN_C_BEGIN
1265e69d4d44SBarry Smith #undef __FUNCT__
1266e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
12677087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1268e69d4d44SBarry Smith {
1269e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1270084e4875SJed Brown   PetscErrorCode  ierr;
1271e69d4d44SBarry Smith 
1272e69d4d44SBarry Smith   PetscFunctionBegin;
1273084e4875SJed Brown   jac->schurpre = ptype;
1274084e4875SJed Brown   if (pre) {
12756bf464f9SBarry Smith     ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
1276084e4875SJed Brown     jac->schur_user = pre;
1277084e4875SJed Brown     ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1278084e4875SJed Brown   }
1279e69d4d44SBarry Smith   PetscFunctionReturn(0);
1280e69d4d44SBarry Smith }
1281e69d4d44SBarry Smith EXTERN_C_END
1282e69d4d44SBarry Smith 
128330ad9308SMatthew Knepley #undef __FUNCT__
128430ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
128530ad9308SMatthew Knepley /*@C
12868c03b21aSDmitry Karpeev    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
128730ad9308SMatthew Knepley 
128830ad9308SMatthew Knepley    Collective on KSP
128930ad9308SMatthew Knepley 
129030ad9308SMatthew Knepley    Input Parameter:
129130ad9308SMatthew Knepley .  pc - the preconditioner context
129230ad9308SMatthew Knepley 
129330ad9308SMatthew Knepley    Output Parameters:
1294a04f6461SBarry Smith +  A00 - the (0,0) block
1295a04f6461SBarry Smith .  A01 - the (0,1) block
1296a04f6461SBarry Smith .  A10 - the (1,0) block
1297a04f6461SBarry Smith -  A11 - the (1,1) block
129830ad9308SMatthew Knepley 
129930ad9308SMatthew Knepley    Level: advanced
130030ad9308SMatthew Knepley 
130130ad9308SMatthew Knepley .seealso: PCFIELDSPLIT
130230ad9308SMatthew Knepley @*/
1303a04f6461SBarry Smith PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
130430ad9308SMatthew Knepley {
130530ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit *) pc->data;
130630ad9308SMatthew Knepley 
130730ad9308SMatthew Knepley   PetscFunctionBegin;
13080700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1309c1235816SBarry Smith   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
1310a04f6461SBarry Smith   if (A00) *A00 = jac->pmat[0];
1311a04f6461SBarry Smith   if (A01) *A01 = jac->B;
1312a04f6461SBarry Smith   if (A10) *A10 = jac->C;
1313a04f6461SBarry Smith   if (A11) *A11 = jac->pmat[1];
131430ad9308SMatthew Knepley   PetscFunctionReturn(0);
131530ad9308SMatthew Knepley }
131630ad9308SMatthew Knepley 
131779416396SBarry Smith EXTERN_C_BEGIN
131879416396SBarry Smith #undef __FUNCT__
131979416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
13207087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
132179416396SBarry Smith {
132279416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1323e69d4d44SBarry Smith   PetscErrorCode ierr;
132479416396SBarry Smith 
132579416396SBarry Smith   PetscFunctionBegin;
132679416396SBarry Smith   jac->type = type;
13273b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
13283b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
13293b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
1330e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1331e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1332e69d4d44SBarry Smith 
13333b224e63SBarry Smith   } else {
13343b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
13353b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
1336e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
13379e7d6b0aSBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr);
13383b224e63SBarry Smith   }
133979416396SBarry Smith   PetscFunctionReturn(0);
134079416396SBarry Smith }
134179416396SBarry Smith EXTERN_C_END
134279416396SBarry Smith 
134351f519a2SBarry Smith EXTERN_C_BEGIN
134451f519a2SBarry Smith #undef __FUNCT__
134551f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
13467087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
134751f519a2SBarry Smith {
134851f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
134951f519a2SBarry Smith 
135051f519a2SBarry Smith   PetscFunctionBegin;
135165e19b50SBarry Smith   if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
135265e19b50SBarry Smith   if (jac->bs > 0 && jac->bs != bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot change fieldsplit blocksize from %D to %D after it has been set",jac->bs,bs);
135351f519a2SBarry Smith   jac->bs = bs;
135451f519a2SBarry Smith   PetscFunctionReturn(0);
135551f519a2SBarry Smith }
135651f519a2SBarry Smith EXTERN_C_END
135751f519a2SBarry Smith 
135879416396SBarry Smith #undef __FUNCT__
135979416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
1360bc08b0f1SBarry Smith /*@
136179416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
136279416396SBarry Smith 
136379416396SBarry Smith    Collective on PC
136479416396SBarry Smith 
136579416396SBarry Smith    Input Parameter:
136679416396SBarry Smith .  pc - the preconditioner context
136781540f2fSBarry Smith .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
136879416396SBarry Smith 
136979416396SBarry Smith    Options Database Key:
1370a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
137179416396SBarry Smith 
137279416396SBarry Smith    Level: Developer
137379416396SBarry Smith 
137479416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
137579416396SBarry Smith 
137679416396SBarry Smith .seealso: PCCompositeSetType()
137779416396SBarry Smith 
137879416396SBarry Smith @*/
13797087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
138079416396SBarry Smith {
13814ac538c5SBarry Smith   PetscErrorCode ierr;
138279416396SBarry Smith 
138379416396SBarry Smith   PetscFunctionBegin;
13840700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13854ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr);
138679416396SBarry Smith   PetscFunctionReturn(0);
138779416396SBarry Smith }
138879416396SBarry Smith 
13890971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
13900971522cSBarry Smith /*MC
1391a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
1392a04f6461SBarry Smith                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.
13930971522cSBarry Smith 
1394edf189efSBarry Smith      To set options on the solvers for each block append -fieldsplit_ to all the PC
1395edf189efSBarry Smith         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
13960971522cSBarry Smith 
1397a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
139869a612a9SBarry Smith          and set the options directly on the resulting KSP object
13990971522cSBarry Smith 
14000971522cSBarry Smith    Level: intermediate
14010971522cSBarry Smith 
140279416396SBarry Smith    Options Database Keys:
140381540f2fSBarry Smith +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
140481540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
140581540f2fSBarry Smith                               been supplied explicitly by -pc_fieldsplit_%d_fields
140681540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
140781540f2fSBarry Smith .   -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative>
1408e69d4d44SBarry Smith .   -pc_fieldsplit_schur_precondition <true,false> default is true
1409435f959eSBarry Smith .   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver
141079416396SBarry Smith 
1411*5d4c12cdSJungho Lee -    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
1412*5d4c12cdSJungho Lee      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
1413*5d4c12cdSJungho Lee 
1414*5d4c12cdSJungho Lee    Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1415d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
1416d32f9abdSBarry Smith 
1417d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
1418d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1419d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
1420d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
1421d32f9abdSBarry Smith 
1422*5d4c12cdSJungho Lee      For the Schur complement preconditioner if J = ( A00 A01 )
1423*5d4c12cdSJungho Lee                                                     ( A10 A11 )
1424*5d4c12cdSJungho Lee      the preconditioner is
1425*5d4c12cdSJungho Lee               (I   -B ksp(A00)) ( inv(A00)   0  (I         0  )
1426*5d4c12cdSJungho Lee               (0    I       ) (   0    ksp(S) ) (-A10 ksp(A00) I  )
1427a04f6461SBarry Smith      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of
1428*5d4c12cdSJungho Lee      ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given in providing the SECOND split or 1 if not give).
1429*5d4c12cdSJungho Lee      For PCFieldSplitGetKSP() when field number is
1430*5d4c12cdSJungho Lee      0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
1431a04f6461SBarry Smith      A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
1432*5d4c12cdSJungho Lee      option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc
1433e69d4d44SBarry Smith 
1434edf189efSBarry Smith      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
1435edf189efSBarry Smith      is used automatically for a second block.
1436edf189efSBarry Smith 
1437ff218e97SBarry Smith      The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1.
1438ff218e97SBarry Smith      Generally it should be used with the AIJ format.
1439ff218e97SBarry Smith 
1440ff218e97SBarry Smith      The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see,
1441ff218e97SBarry Smith      for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT
1442ff218e97SBarry Smith      inside a smoother resulting in "Distributive Smoothers".
14430716a85fSBarry Smith 
1444a541d17aSBarry Smith    Concepts: physics based preconditioners, block preconditioners
14450971522cSBarry Smith 
14467e8cb189SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
1447e69d4d44SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
14480971522cSBarry Smith M*/
14490971522cSBarry Smith 
14500971522cSBarry Smith EXTERN_C_BEGIN
14510971522cSBarry Smith #undef __FUNCT__
14520971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
14537087cfbeSBarry Smith PetscErrorCode  PCCreate_FieldSplit(PC pc)
14540971522cSBarry Smith {
14550971522cSBarry Smith   PetscErrorCode ierr;
14560971522cSBarry Smith   PC_FieldSplit  *jac;
14570971522cSBarry Smith 
14580971522cSBarry Smith   PetscFunctionBegin;
145938f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
14600971522cSBarry Smith   jac->bs        = -1;
14610971522cSBarry Smith   jac->nsplits   = 0;
14623e197d65SBarry Smith   jac->type      = PC_COMPOSITE_MULTIPLICATIVE;
1463e6cab6aaSJed Brown   jac->schurpre  = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
1464c5d2311dSJed Brown   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL;
146551f519a2SBarry Smith 
14660971522cSBarry Smith   pc->data     = (void*)jac;
14670971522cSBarry Smith 
14680971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
1469421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
14700971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
1471574deadeSBarry Smith   pc->ops->reset             = PCReset_FieldSplit;
14720971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
14730971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
14740971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
14750971522cSBarry Smith   pc->ops->applyrichardson   = 0;
14760971522cSBarry Smith 
147769a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
147869a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
14790971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
14800971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1481704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
1482704ba839SBarry Smith                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
148379416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
148479416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
148551f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
148651f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
14870971522cSBarry Smith   PetscFunctionReturn(0);
14880971522cSBarry Smith }
14890971522cSBarry Smith EXTERN_C_END
14900971522cSBarry Smith 
14910971522cSBarry Smith 
1492a541d17aSBarry Smith 
1493