xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision ad4df1005a67cf9575a5c4294f11fd96229b2bb4)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
36356e834SBarry Smith #include "private/pcimpl.h"     /*I "petscpc.h" I*/
40971522cSBarry Smith 
5c5d2311dSJed Brown const char *PCFieldSplitSchurPreTypes[] = {"SELF","DIAG","USER","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0};
6c5d2311dSJed Brown const char *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;
210971522cSBarry Smith   PetscInt          *fields;
221b9fc7fcSBarry Smith   VecScatter        sctx;
234aa3045dSJed Brown   IS                is;
2451f519a2SBarry Smith   PC_FieldSplitLink next,previous;
250971522cSBarry Smith };
260971522cSBarry Smith 
270971522cSBarry Smith typedef struct {
2868dd23aaSBarry Smith   PCCompositeType   type;
29a4efd8eaSMatthew Knepley   PetscTruth        defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */
306c924f48SJed Brown   PetscTruth        splitdefined; /* Flag is set after the splits have been defined, to prevent more splits from being added */
31519d70e2SJed Brown   PetscTruth        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 */
38704ba839SBarry Smith   PetscTruth        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 */
4230ad9308SMatthew Knepley   Mat               schur;        /* The Schur complement S = D - C A^{-1} B */
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;
480971522cSBarry Smith } PC_FieldSplit;
490971522cSBarry Smith 
5016913363SBarry Smith /*
5116913363SBarry Smith     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
5216913363SBarry Smith    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
5316913363SBarry Smith    PC you could change this.
5416913363SBarry Smith */
55084e4875SJed Brown 
56e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it.  This way the
57084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */
58084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac)
59084e4875SJed Brown {
60084e4875SJed Brown   switch (jac->schurpre) {
61084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur;
62084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1];
63084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */
64084e4875SJed Brown     default:
65084e4875SJed Brown       return jac->schur_user ? jac->schur_user : jac->pmat[1];
66084e4875SJed Brown   }
67084e4875SJed Brown }
68084e4875SJed Brown 
69084e4875SJed Brown 
700971522cSBarry Smith #undef __FUNCT__
710971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
720971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
730971522cSBarry Smith {
740971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
750971522cSBarry Smith   PetscErrorCode    ierr;
760971522cSBarry Smith   PetscTruth        iascii;
770971522cSBarry Smith   PetscInt          i,j;
785a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
790971522cSBarry Smith 
800971522cSBarry Smith   PetscFunctionBegin;
812692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
820971522cSBarry Smith   if (iascii) {
8351f519a2SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
8469a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
850971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
860971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
871ab39975SBarry Smith       if (ilink->fields) {
880971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
8979416396SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
905a9f2f41SSatish Balay 	for (j=0; j<ilink->nfields; j++) {
9179416396SBarry Smith 	  if (j > 0) {
9279416396SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
9379416396SBarry Smith 	  }
945a9f2f41SSatish Balay 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
950971522cSBarry Smith 	}
960971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
9779416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
981ab39975SBarry Smith       } else {
991ab39975SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1001ab39975SBarry Smith       }
1015a9f2f41SSatish Balay       ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
1025a9f2f41SSatish Balay       ilink = ilink->next;
1030971522cSBarry Smith     }
1040971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1050971522cSBarry Smith   } else {
10665e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1070971522cSBarry Smith   }
1080971522cSBarry Smith   PetscFunctionReturn(0);
1090971522cSBarry Smith }
1100971522cSBarry Smith 
1110971522cSBarry Smith #undef __FUNCT__
1123b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur"
1133b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
1143b224e63SBarry Smith {
1153b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1163b224e63SBarry Smith   PetscErrorCode    ierr;
1173b224e63SBarry Smith   PetscTruth        iascii;
1183b224e63SBarry Smith   PetscInt          i,j;
1193b224e63SBarry Smith   PC_FieldSplitLink ilink = jac->head;
1203b224e63SBarry Smith   KSP               ksp;
1213b224e63SBarry Smith 
1223b224e63SBarry Smith   PetscFunctionBegin;
1232692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1243b224e63SBarry Smith   if (iascii) {
125c5d2311dSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactorizationTypes[jac->schurfactorization]);CHKERRQ(ierr);
1263b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Split info:\n");CHKERRQ(ierr);
1273b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1283b224e63SBarry Smith     for (i=0; i<jac->nsplits; i++) {
1293b224e63SBarry Smith       if (ilink->fields) {
1303b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
1313b224e63SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
1323b224e63SBarry Smith 	for (j=0; j<ilink->nfields; j++) {
1333b224e63SBarry Smith 	  if (j > 0) {
1343b224e63SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
1353b224e63SBarry Smith 	  }
1363b224e63SBarry Smith 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
1373b224e63SBarry Smith 	}
1383b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1393b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1403b224e63SBarry Smith       } else {
1413b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1423b224e63SBarry Smith       }
1433b224e63SBarry Smith       ilink = ilink->next;
1443b224e63SBarry Smith     }
1453b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A block \n");CHKERRQ(ierr);
1463b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
14712cae6f2SJed Brown     if (jac->schur) {
14812cae6f2SJed Brown       ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
1493b224e63SBarry Smith       ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
15012cae6f2SJed Brown     } else {
15112cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
15212cae6f2SJed Brown     }
1533b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1543b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = D - C inv(A) B \n");CHKERRQ(ierr);
1553b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15612cae6f2SJed Brown     if (jac->kspschur) {
1573b224e63SBarry Smith       ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
15812cae6f2SJed Brown     } else {
15912cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
16012cae6f2SJed Brown     }
1613b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1623b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1633b224e63SBarry Smith   } else {
16465e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1653b224e63SBarry Smith   }
1663b224e63SBarry Smith   PetscFunctionReturn(0);
1673b224e63SBarry Smith }
1683b224e63SBarry Smith 
1693b224e63SBarry Smith #undef __FUNCT__
1706c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private"
1716c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */
1726c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc)
1736c924f48SJed Brown {
1746c924f48SJed Brown   PetscErrorCode ierr;
1756c924f48SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1766c924f48SJed Brown   PetscInt       i,nfields,*ifields;
1776c924f48SJed Brown   PetscTruth     flg;
1786c924f48SJed Brown   char           optionname[128],splitname[8];
1796c924f48SJed Brown 
1806c924f48SJed Brown   PetscFunctionBegin;
1816c924f48SJed Brown   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr);
1826c924f48SJed Brown   for (i=0,flg=PETSC_TRUE; ; i++) {
1836c924f48SJed Brown     ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
1846c924f48SJed Brown     ierr = PetscSNPrintf(optionname,sizeof optionname,"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr);
1856c924f48SJed Brown     nfields = jac->bs;
1866c924f48SJed Brown     ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr);
1876c924f48SJed Brown     if (!flg) break;
1886c924f48SJed Brown     if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
1896c924f48SJed Brown     ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields);CHKERRQ(ierr);
1906c924f48SJed Brown   }
1916c924f48SJed Brown   if (i > 0) {
1926c924f48SJed Brown     /* Makes command-line setting of splits take precedence over setting them in code.
1936c924f48SJed Brown        Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would
1946c924f48SJed Brown        create new splits, which would probably not be what the user wanted. */
1956c924f48SJed Brown     jac->splitdefined = PETSC_TRUE;
1966c924f48SJed Brown   }
1976c924f48SJed Brown   ierr = PetscFree(ifields);CHKERRQ(ierr);
1986c924f48SJed Brown   PetscFunctionReturn(0);
1996c924f48SJed Brown }
2006c924f48SJed Brown 
2016c924f48SJed Brown #undef __FUNCT__
20269a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
20369a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
2040971522cSBarry Smith {
2050971522cSBarry Smith   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
2060971522cSBarry Smith   PetscErrorCode    ierr;
2075a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
2086c924f48SJed Brown   PetscTruth        flg = PETSC_FALSE;
2096c924f48SJed Brown   PetscInt          i;
2100971522cSBarry Smith 
2110971522cSBarry Smith   PetscFunctionBegin;
212d32f9abdSBarry Smith   if (!ilink) {
213d32f9abdSBarry Smith 
214521d7252SBarry Smith     if (jac->bs <= 0) {
215704ba839SBarry Smith       if (pc->pmat) {
216521d7252SBarry Smith         ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
217704ba839SBarry Smith       } else {
218704ba839SBarry Smith         jac->bs = 1;
219704ba839SBarry Smith       }
220521d7252SBarry Smith     }
221d32f9abdSBarry Smith 
222d32f9abdSBarry Smith     ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
223d32f9abdSBarry Smith     if (!flg) {
224d32f9abdSBarry Smith       /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
225d32f9abdSBarry Smith          then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
2266c924f48SJed Brown       ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
2276c924f48SJed Brown       if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
228d32f9abdSBarry Smith     }
2296c924f48SJed Brown     if (flg || !jac->splitdefined) {
230d32f9abdSBarry Smith       ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
231db4c96c1SJed Brown       for (i=0; i<jac->bs; i++) {
2326c924f48SJed Brown         char splitname[8];
2336c924f48SJed Brown         ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
234db4c96c1SJed Brown         ierr = PCFieldSplitSetFields(pc,splitname,1,&i);CHKERRQ(ierr);
23579416396SBarry Smith       }
23697bbdb24SBarry Smith       jac->defaultsplit = PETSC_TRUE;
237521d7252SBarry Smith     }
238edf189efSBarry Smith   } else if (jac->nsplits == 1) {
239edf189efSBarry Smith     if (ilink->is) {
240edf189efSBarry Smith       IS       is2;
241edf189efSBarry Smith       PetscInt nmin,nmax;
242edf189efSBarry Smith 
243edf189efSBarry Smith       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
244edf189efSBarry Smith       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
245db4c96c1SJed Brown       ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr);
246edf189efSBarry Smith       ierr = ISDestroy(is2);CHKERRQ(ierr);
247e7e72b3dSBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
248edf189efSBarry Smith   }
249e7e72b3dSBarry Smith   if (jac->nsplits < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields");
25069a612a9SBarry Smith   PetscFunctionReturn(0);
25169a612a9SBarry Smith }
25269a612a9SBarry Smith 
25369a612a9SBarry Smith 
25469a612a9SBarry Smith #undef __FUNCT__
25569a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
25669a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
25769a612a9SBarry Smith {
25869a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
25969a612a9SBarry Smith   PetscErrorCode    ierr;
2605a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
261cf502942SBarry Smith   PetscInt          i,nsplit,ccsize;
26269a612a9SBarry Smith   MatStructure      flag = pc->flag;
2636c8605c2SJed Brown   PetscTruth        sorted;
26469a612a9SBarry Smith 
26569a612a9SBarry Smith   PetscFunctionBegin;
26669a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
26797bbdb24SBarry Smith   nsplit = jac->nsplits;
2685a9f2f41SSatish Balay   ilink  = jac->head;
26997bbdb24SBarry Smith 
27097bbdb24SBarry Smith   /* get the matrices for each split */
271704ba839SBarry Smith   if (!jac->issetup) {
2721b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
27397bbdb24SBarry Smith 
274704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
275704ba839SBarry Smith 
276704ba839SBarry Smith     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
27751f519a2SBarry Smith     bs     = jac->bs;
27897bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
279cf502942SBarry Smith     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
2801b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
2811b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
2821b9fc7fcSBarry Smith       if (jac->defaultsplit) {
283704ba839SBarry Smith         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
284704ba839SBarry Smith       } else if (!ilink->is) {
285ccb205f8SBarry Smith         if (ilink->nfields > 1) {
2865a9f2f41SSatish Balay           PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
2875a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
2881b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
2891b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
2901b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
29197bbdb24SBarry Smith             }
29297bbdb24SBarry Smith           }
293704ba839SBarry Smith           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr);
294ccb205f8SBarry Smith           ierr = PetscFree(ii);CHKERRQ(ierr);
295ccb205f8SBarry Smith         } else {
296704ba839SBarry Smith           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
297ccb205f8SBarry Smith         }
2983e197d65SBarry Smith       }
299704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
300e32f2f54SBarry Smith       if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
301704ba839SBarry Smith       ilink = ilink->next;
3021b9fc7fcSBarry Smith     }
3031b9fc7fcSBarry Smith   }
3041b9fc7fcSBarry Smith 
305704ba839SBarry Smith   ilink  = jac->head;
30697bbdb24SBarry Smith   if (!jac->pmat) {
307cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
308cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
3094aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
310704ba839SBarry Smith       ilink = ilink->next;
311cf502942SBarry Smith     }
31297bbdb24SBarry Smith   } else {
313cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
3144aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
315704ba839SBarry Smith       ilink = ilink->next;
316cf502942SBarry Smith     }
31797bbdb24SBarry Smith   }
318519d70e2SJed Brown   if (jac->realdiagonal) {
319519d70e2SJed Brown     ilink = jac->head;
320519d70e2SJed Brown     if (!jac->mat) {
321519d70e2SJed Brown       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr);
322519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
323519d70e2SJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
324519d70e2SJed Brown         ilink = ilink->next;
325519d70e2SJed Brown       }
326519d70e2SJed Brown     } else {
327519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
328519d70e2SJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
329519d70e2SJed Brown         ilink = ilink->next;
330519d70e2SJed Brown       }
331519d70e2SJed Brown     }
332519d70e2SJed Brown   } else {
333519d70e2SJed Brown     jac->mat = jac->pmat;
334519d70e2SJed Brown   }
33597bbdb24SBarry Smith 
3366c8605c2SJed Brown   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
33768dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
33868dd23aaSBarry Smith     ilink  = jac->head;
33968dd23aaSBarry Smith     if (!jac->Afield) {
34068dd23aaSBarry Smith       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
34168dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
3424aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
34368dd23aaSBarry Smith         ilink = ilink->next;
34468dd23aaSBarry Smith       }
34568dd23aaSBarry Smith     } else {
34668dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
3474aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
34868dd23aaSBarry Smith         ilink = ilink->next;
34968dd23aaSBarry Smith       }
35068dd23aaSBarry Smith     }
35168dd23aaSBarry Smith   }
35268dd23aaSBarry Smith 
3533b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
3543b224e63SBarry Smith     IS       ccis;
3554aa3045dSJed Brown     PetscInt rstart,rend;
356e7e72b3dSBarry Smith     if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
35768dd23aaSBarry Smith 
358e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
359e6cab6aaSJed Brown     ierr  = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
360e6cab6aaSJed Brown 
3613b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
3623b224e63SBarry Smith     if (jac->schur) {
3633b224e63SBarry Smith       ilink = jac->head;
3644aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
3654aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
3663b224e63SBarry Smith       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
3673b224e63SBarry Smith       ilink = ilink->next;
3684aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
3694aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
3703b224e63SBarry Smith       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
371519d70e2SJed Brown       ierr  = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr);
372084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr);
3733b224e63SBarry Smith 
3743b224e63SBarry Smith      } else {
3751cee3971SBarry Smith       KSP ksp;
3766c924f48SJed Brown       char schurprefix[256];
3773b224e63SBarry Smith 
3783b224e63SBarry Smith       /* extract the B and C matrices */
3793b224e63SBarry Smith       ilink = jac->head;
3804aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
3814aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
3823b224e63SBarry Smith       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
3833b224e63SBarry Smith       ilink = ilink->next;
3844aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
3854aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
3863b224e63SBarry Smith       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
387084e4875SJed Brown       /* Better would be to use 'mat[0]' (diagonal block of the real matrix) preconditioned by pmat[0] */
388519d70e2SJed Brown       ierr  = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr);
3891cee3971SBarry Smith       ierr  = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
390e69d4d44SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr);
3913b224e63SBarry Smith       ierr  = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
3923b224e63SBarry Smith 
3933b224e63SBarry Smith       ierr  = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr);
3941cee3971SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
395084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
396084e4875SJed Brown       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
397084e4875SJed Brown         PC pc;
398cd5f4a64SJed Brown         ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr);
399084e4875SJed Brown         ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
400084e4875SJed Brown         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
401e69d4d44SBarry Smith       }
4026c924f48SJed Brown       ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
4036c924f48SJed Brown       ierr  = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr);
4043b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
4053b224e63SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
4063b224e63SBarry Smith 
4073b224e63SBarry Smith       ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr);
4083b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr);
4093b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr);
4103b224e63SBarry Smith       ilink = jac->head;
4113b224e63SBarry Smith       ilink->x = jac->x[0]; ilink->y = jac->y[0];
4123b224e63SBarry Smith       ilink = ilink->next;
4133b224e63SBarry Smith       ilink->x = jac->x[1]; ilink->y = jac->y[1];
4143b224e63SBarry Smith     }
4153b224e63SBarry Smith   } else {
41697bbdb24SBarry Smith     /* set up the individual PCs */
41797bbdb24SBarry Smith     i    = 0;
4185a9f2f41SSatish Balay     ilink = jac->head;
4195a9f2f41SSatish Balay     while (ilink) {
420519d70e2SJed Brown       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr);
4213b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
4225a9f2f41SSatish Balay       ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);
4235a9f2f41SSatish Balay       ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
42497bbdb24SBarry Smith       i++;
4255a9f2f41SSatish Balay       ilink = ilink->next;
4260971522cSBarry Smith     }
42797bbdb24SBarry Smith 
42897bbdb24SBarry Smith     /* create work vectors for each split */
4291b9fc7fcSBarry Smith     if (!jac->x) {
43097bbdb24SBarry Smith       ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
4315a9f2f41SSatish Balay       ilink = jac->head;
43297bbdb24SBarry Smith       for (i=0; i<nsplit; i++) {
433906ed7ccSBarry Smith         Vec *vl,*vr;
434906ed7ccSBarry Smith 
435906ed7ccSBarry Smith         ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
436906ed7ccSBarry Smith         ilink->x  = *vr;
437906ed7ccSBarry Smith         ilink->y  = *vl;
438906ed7ccSBarry Smith         ierr      = PetscFree(vr);CHKERRQ(ierr);
439906ed7ccSBarry Smith         ierr      = PetscFree(vl);CHKERRQ(ierr);
4405a9f2f41SSatish Balay         jac->x[i] = ilink->x;
4415a9f2f41SSatish Balay         jac->y[i] = ilink->y;
4425a9f2f41SSatish Balay         ilink     = ilink->next;
44397bbdb24SBarry Smith       }
4443b224e63SBarry Smith     }
4453b224e63SBarry Smith   }
4463b224e63SBarry Smith 
4473b224e63SBarry Smith 
448d0f46423SBarry Smith   if (!jac->head->sctx) {
4493b224e63SBarry Smith     Vec xtmp;
4503b224e63SBarry Smith 
45179416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
4521b9fc7fcSBarry Smith 
4535a9f2f41SSatish Balay     ilink = jac->head;
4541b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
4551b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
456704ba839SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
4575a9f2f41SSatish Balay       ilink = ilink->next;
4581b9fc7fcSBarry Smith     }
4591b9fc7fcSBarry Smith     ierr = VecDestroy(xtmp);CHKERRQ(ierr);
4601b9fc7fcSBarry Smith   }
4610971522cSBarry Smith   PetscFunctionReturn(0);
4620971522cSBarry Smith }
4630971522cSBarry Smith 
4645a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
465ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
466ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
4675a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
468ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
469ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
47079416396SBarry Smith 
4710971522cSBarry Smith #undef __FUNCT__
4723b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur"
4733b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
4743b224e63SBarry Smith {
4753b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
4763b224e63SBarry Smith   PetscErrorCode    ierr;
4773b224e63SBarry Smith   KSP               ksp;
4783b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
4793b224e63SBarry Smith 
4803b224e63SBarry Smith   PetscFunctionBegin;
4813b224e63SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
4823b224e63SBarry Smith 
483c5d2311dSJed Brown   switch (jac->schurfactorization) {
484c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG:
485c5d2311dSJed Brown       /* [A 0; 0 -S], positive definite, suitable for MINRES */
486c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
487c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
488c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
489c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
490c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
491c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
492c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
493c5d2311dSJed Brown       ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr);
494c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
495c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
496c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
497c5d2311dSJed Brown       break;
498c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER:
499c5d2311dSJed Brown       /* [A 0; C S], suitable for left preconditioning */
500c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
501c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
502c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
503c5d2311dSJed Brown       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
504c5d2311dSJed Brown       ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr);
505c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
506c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
507c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
508c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
509c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
510c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
511c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
512c5d2311dSJed Brown       break;
513c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER:
514c5d2311dSJed Brown       /* [A B; 0 S], suitable for right preconditioning */
515c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
516c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
517c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
518c5d2311dSJed Brown       ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
519c5d2311dSJed Brown       ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr);
520c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
521c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
522c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
523c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
524c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
525c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
526c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
527c5d2311dSJed Brown       break;
528c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL:
529c5d2311dSJed Brown       /* [1 0; CA^{-1} 1] [A 0; 0 S] [1 A^{-1}B; 0 1], an exact solve if applied exactly, needs one extra solve with A */
5303b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5313b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5323b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
5333b224e63SBarry Smith       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
5343b224e63SBarry Smith       ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
5353b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5363b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5373b224e63SBarry Smith 
5383b224e63SBarry Smith       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
5393b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5403b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5413b224e63SBarry Smith 
5423b224e63SBarry Smith       ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
5433b224e63SBarry Smith       ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
5443b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
5453b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5463b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
547c5d2311dSJed Brown   }
5483b224e63SBarry Smith   PetscFunctionReturn(0);
5493b224e63SBarry Smith }
5503b224e63SBarry Smith 
5513b224e63SBarry Smith #undef __FUNCT__
5520971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
5530971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
5540971522cSBarry Smith {
5550971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
5560971522cSBarry Smith   PetscErrorCode    ierr;
5575a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
558af93645bSJed Brown   PetscInt          cnt;
5590971522cSBarry Smith 
5600971522cSBarry Smith   PetscFunctionBegin;
56151f519a2SBarry Smith   CHKMEMQ;
56251f519a2SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
56351f519a2SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
56451f519a2SBarry Smith 
56579416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
5661b9fc7fcSBarry Smith     if (jac->defaultsplit) {
5670971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
5685a9f2f41SSatish Balay       while (ilink) {
5695a9f2f41SSatish Balay         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
5705a9f2f41SSatish Balay         ilink = ilink->next;
5710971522cSBarry Smith       }
5720971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
5731b9fc7fcSBarry Smith     } else {
574efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
5755a9f2f41SSatish Balay       while (ilink) {
5765a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
5775a9f2f41SSatish Balay         ilink = ilink->next;
5781b9fc7fcSBarry Smith       }
5791b9fc7fcSBarry Smith     }
58016913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
58179416396SBarry Smith     if (!jac->w1) {
58279416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
58379416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
58479416396SBarry Smith     }
585efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
5865a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
5873e197d65SBarry Smith     cnt = 1;
5885a9f2f41SSatish Balay     while (ilink->next) {
5895a9f2f41SSatish Balay       ilink = ilink->next;
5903e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
5913e197d65SBarry Smith       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
5923e197d65SBarry Smith       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
5933e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5943e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5953e197d65SBarry Smith       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
5963e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5973e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5983e197d65SBarry Smith     }
59951f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
60011755939SBarry Smith       cnt -= 2;
60151f519a2SBarry Smith       while (ilink->previous) {
60251f519a2SBarry Smith         ilink = ilink->previous;
60311755939SBarry Smith         /* compute the residual only over the part of the vector needed */
60411755939SBarry Smith         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
60511755939SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
60611755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
60711755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
60811755939SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
60911755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
61011755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
61151f519a2SBarry Smith       }
61211755939SBarry Smith     }
61365e19b50SBarry Smith   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
61451f519a2SBarry Smith   CHKMEMQ;
6150971522cSBarry Smith   PetscFunctionReturn(0);
6160971522cSBarry Smith }
6170971522cSBarry Smith 
618421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
619ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
620ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
621421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
622ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
623ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
624421e10b8SBarry Smith 
625421e10b8SBarry Smith #undef __FUNCT__
626421e10b8SBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
627421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
628421e10b8SBarry Smith {
629421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
630421e10b8SBarry Smith   PetscErrorCode    ierr;
631421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
632421e10b8SBarry Smith 
633421e10b8SBarry Smith   PetscFunctionBegin;
634421e10b8SBarry Smith   CHKMEMQ;
635421e10b8SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
636421e10b8SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
637421e10b8SBarry Smith 
638421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
639421e10b8SBarry Smith     if (jac->defaultsplit) {
640421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
641421e10b8SBarry Smith       while (ilink) {
642421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
643421e10b8SBarry Smith 	ilink = ilink->next;
644421e10b8SBarry Smith       }
645421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
646421e10b8SBarry Smith     } else {
647421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
648421e10b8SBarry Smith       while (ilink) {
649421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
650421e10b8SBarry Smith 	ilink = ilink->next;
651421e10b8SBarry Smith       }
652421e10b8SBarry Smith     }
653421e10b8SBarry Smith   } else {
654421e10b8SBarry Smith     if (!jac->w1) {
655421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
656421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
657421e10b8SBarry Smith     }
658421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
659421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
660421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
661421e10b8SBarry Smith       while (ilink->next) {
662421e10b8SBarry Smith         ilink = ilink->next;
6639989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
664421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
665421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
666421e10b8SBarry Smith       }
667421e10b8SBarry Smith       while (ilink->previous) {
668421e10b8SBarry Smith         ilink = ilink->previous;
6699989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
670421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
671421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
672421e10b8SBarry Smith       }
673421e10b8SBarry Smith     } else {
674421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
675421e10b8SBarry Smith 	ilink = ilink->next;
676421e10b8SBarry Smith       }
677421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
678421e10b8SBarry Smith       while (ilink->previous) {
679421e10b8SBarry Smith 	ilink = ilink->previous;
6809989ab13SBarry Smith 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
681421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
682421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
683421e10b8SBarry Smith       }
684421e10b8SBarry Smith     }
685421e10b8SBarry Smith   }
686421e10b8SBarry Smith   CHKMEMQ;
687421e10b8SBarry Smith   PetscFunctionReturn(0);
688421e10b8SBarry Smith }
689421e10b8SBarry Smith 
6900971522cSBarry Smith #undef __FUNCT__
6910971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
6920971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
6930971522cSBarry Smith {
6940971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
6950971522cSBarry Smith   PetscErrorCode    ierr;
6965a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
6970971522cSBarry Smith 
6980971522cSBarry Smith   PetscFunctionBegin;
6995a9f2f41SSatish Balay   while (ilink) {
7005a9f2f41SSatish Balay     ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr);
7015a9f2f41SSatish Balay     if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);}
7025a9f2f41SSatish Balay     if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);}
7035a9f2f41SSatish Balay     if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);}
704704ba839SBarry Smith     if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);}
7055a9f2f41SSatish Balay     next = ilink->next;
7067e8c30b6SJed Brown     ierr = PetscFree(ilink->splitname);CHKERRQ(ierr);
707704ba839SBarry Smith     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
708704ba839SBarry Smith     ierr = PetscFree(ilink);CHKERRQ(ierr);
7095a9f2f41SSatish Balay     ilink = next;
7100971522cSBarry Smith   }
71105b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
712519d70e2SJed Brown   if (jac->mat && jac->mat != jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);}
71397bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
71468dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
71579416396SBarry Smith   if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);}
71679416396SBarry Smith   if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);}
7173b224e63SBarry Smith   if (jac->schur) {ierr = MatDestroy(jac->schur);CHKERRQ(ierr);}
718084e4875SJed Brown   if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);}
719d0f46423SBarry Smith   if (jac->kspschur) {ierr = KSPDestroy(jac->kspschur);CHKERRQ(ierr);}
7203b224e63SBarry Smith   if (jac->B) {ierr = MatDestroy(jac->B);CHKERRQ(ierr);}
7213b224e63SBarry Smith   if (jac->C) {ierr = MatDestroy(jac->C);CHKERRQ(ierr);}
7220971522cSBarry Smith   ierr = PetscFree(jac);CHKERRQ(ierr);
7230971522cSBarry Smith   PetscFunctionReturn(0);
7240971522cSBarry Smith }
7250971522cSBarry Smith 
7260971522cSBarry Smith #undef __FUNCT__
7270971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
7280971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
7290971522cSBarry Smith {
7301b9fc7fcSBarry Smith   PetscErrorCode  ierr;
7316c924f48SJed Brown   PetscInt        bs;
732084e4875SJed Brown   PetscTruth      flg;
7339dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
7343b224e63SBarry Smith   PCCompositeType ctype;
7351b9fc7fcSBarry Smith 
7360971522cSBarry Smith   PetscFunctionBegin;
7371b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
738519d70e2SJed Brown   ierr = PetscOptionsTruth("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr);
73951f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
74051f519a2SBarry Smith   if (flg) {
74151f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
74251f519a2SBarry Smith   }
743704ba839SBarry Smith 
7443b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
7453b224e63SBarry Smith   if (flg) {
7463b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
7473b224e63SBarry Smith   }
748d32f9abdSBarry Smith 
749c30613efSMatthew Knepley   /* Only setup fields once */
750c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
751d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
752d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
7536c924f48SJed Brown     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
7546c924f48SJed Brown     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
755d32f9abdSBarry Smith   }
756c5d2311dSJed Brown   if (jac->type == PC_COMPOSITE_SCHUR) {
757c5d2311dSJed Brown     ierr = PetscOptionsEnum("-pc_fieldsplit_schur_factorization_type","Factorization to use","None",PCFieldSplitSchurFactorizationTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,PETSC_NULL);CHKERRQ(ierr);
758084e4875SJed 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);
759c5d2311dSJed Brown   }
7601b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
7610971522cSBarry Smith   PetscFunctionReturn(0);
7620971522cSBarry Smith }
7630971522cSBarry Smith 
7640971522cSBarry Smith /*------------------------------------------------------------------------------------*/
7650971522cSBarry Smith 
7660971522cSBarry Smith EXTERN_C_BEGIN
7670971522cSBarry Smith #undef __FUNCT__
7680971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
7696685144eSJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields)
7700971522cSBarry Smith {
77197bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
7720971522cSBarry Smith   PetscErrorCode    ierr;
7735a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
77469a612a9SBarry Smith   char              prefix[128];
77551f519a2SBarry Smith   PetscInt          i;
7760971522cSBarry Smith 
7770971522cSBarry Smith   PetscFunctionBegin;
7786c924f48SJed Brown   if (jac->splitdefined) {
7796c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
7806c924f48SJed Brown     PetscFunctionReturn(0);
7816c924f48SJed Brown   }
78251f519a2SBarry Smith   for (i=0; i<n; i++) {
783e32f2f54SBarry 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);
784e32f2f54SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
78551f519a2SBarry Smith   }
786704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
787db4c96c1SJed Brown   ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
788704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
7895a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
7905a9f2f41SSatish Balay   ilink->nfields = n;
7915a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
7927adad957SLisandro Dalcin   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
7931cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
7945a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
79569a612a9SBarry Smith 
7966c924f48SJed Brown   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",splitname);CHKERRQ(ierr);
7975a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
7980971522cSBarry Smith 
7990971522cSBarry Smith   if (!next) {
8005a9f2f41SSatish Balay     jac->head       = ilink;
80151f519a2SBarry Smith     ilink->previous = PETSC_NULL;
8020971522cSBarry Smith   } else {
8030971522cSBarry Smith     while (next->next) {
8040971522cSBarry Smith       next = next->next;
8050971522cSBarry Smith     }
8065a9f2f41SSatish Balay     next->next      = ilink;
80751f519a2SBarry Smith     ilink->previous = next;
8080971522cSBarry Smith   }
8090971522cSBarry Smith   jac->nsplits++;
8100971522cSBarry Smith   PetscFunctionReturn(0);
8110971522cSBarry Smith }
8120971522cSBarry Smith EXTERN_C_END
8130971522cSBarry Smith 
814e69d4d44SBarry Smith EXTERN_C_BEGIN
815e69d4d44SBarry Smith #undef __FUNCT__
816e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
817e69d4d44SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
818e69d4d44SBarry Smith {
819e69d4d44SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
820e69d4d44SBarry Smith   PetscErrorCode ierr;
821e69d4d44SBarry Smith 
822e69d4d44SBarry Smith   PetscFunctionBegin;
823519d70e2SJed Brown   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
824e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
825e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
826084e4875SJed Brown   *n = jac->nsplits;
827e69d4d44SBarry Smith   PetscFunctionReturn(0);
828e69d4d44SBarry Smith }
829e69d4d44SBarry Smith EXTERN_C_END
8300971522cSBarry Smith 
8310971522cSBarry Smith EXTERN_C_BEGIN
8320971522cSBarry Smith #undef __FUNCT__
83369a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
834dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
8350971522cSBarry Smith {
8360971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
8370971522cSBarry Smith   PetscErrorCode    ierr;
8380971522cSBarry Smith   PetscInt          cnt = 0;
8395a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
8400971522cSBarry Smith 
8410971522cSBarry Smith   PetscFunctionBegin;
84269a612a9SBarry Smith   ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr);
8435a9f2f41SSatish Balay   while (ilink) {
8445a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
8455a9f2f41SSatish Balay     ilink = ilink->next;
8460971522cSBarry Smith   }
847e32f2f54SBarry Smith   if (cnt != jac->nsplits) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits);
8480971522cSBarry Smith   *n = jac->nsplits;
8490971522cSBarry Smith   PetscFunctionReturn(0);
8500971522cSBarry Smith }
8510971522cSBarry Smith EXTERN_C_END
8520971522cSBarry Smith 
853704ba839SBarry Smith EXTERN_C_BEGIN
854704ba839SBarry Smith #undef __FUNCT__
855704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
856db4c96c1SJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
857704ba839SBarry Smith {
858704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
859704ba839SBarry Smith   PetscErrorCode    ierr;
860704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
861704ba839SBarry Smith   char              prefix[128];
862704ba839SBarry Smith 
863704ba839SBarry Smith   PetscFunctionBegin;
8646c924f48SJed Brown   if (jac->splitdefined) {
8656c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
8666c924f48SJed Brown     PetscFunctionReturn(0);
8676c924f48SJed Brown   }
86816913363SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
869db4c96c1SJed Brown   ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
8701ab39975SBarry Smith   ilink->is      = is;
8711ab39975SBarry Smith   ierr           = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
872704ba839SBarry Smith   ilink->next    = PETSC_NULL;
873704ba839SBarry Smith   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
8741cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
875704ba839SBarry Smith   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
876704ba839SBarry Smith 
8776c924f48SJed Brown   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",splitname);CHKERRQ(ierr);
878704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
879704ba839SBarry Smith 
880704ba839SBarry Smith   if (!next) {
881704ba839SBarry Smith     jac->head       = ilink;
882704ba839SBarry Smith     ilink->previous = PETSC_NULL;
883704ba839SBarry Smith   } else {
884704ba839SBarry Smith     while (next->next) {
885704ba839SBarry Smith       next = next->next;
886704ba839SBarry Smith     }
887704ba839SBarry Smith     next->next      = ilink;
888704ba839SBarry Smith     ilink->previous = next;
889704ba839SBarry Smith   }
890704ba839SBarry Smith   jac->nsplits++;
891704ba839SBarry Smith 
892704ba839SBarry Smith   PetscFunctionReturn(0);
893704ba839SBarry Smith }
894704ba839SBarry Smith EXTERN_C_END
895704ba839SBarry Smith 
8960971522cSBarry Smith #undef __FUNCT__
8970971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
8980971522cSBarry Smith /*@
8990971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
9000971522cSBarry Smith 
901*ad4df100SBarry Smith     Logically Collective on PC
9020971522cSBarry Smith 
9030971522cSBarry Smith     Input Parameters:
9040971522cSBarry Smith +   pc  - the preconditioner context
905db4c96c1SJed Brown .   splitname - name of this split
9060971522cSBarry Smith .   n - the number of fields in this split
907db4c96c1SJed Brown -   fields - the fields in this split
9080971522cSBarry Smith 
9090971522cSBarry Smith     Level: intermediate
9100971522cSBarry Smith 
911d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
912d32f9abdSBarry Smith 
913d32f9abdSBarry Smith      The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block
914d32f9abdSBarry 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
915d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
916d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
917d32f9abdSBarry Smith 
918db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
919db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
920db4c96c1SJed Brown 
921d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
9220971522cSBarry Smith 
9230971522cSBarry Smith @*/
9246685144eSJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields)
9250971522cSBarry Smith {
9266685144eSJed Brown   PetscErrorCode ierr,(*f)(PC,const char[],PetscInt,const PetscInt *);
9270971522cSBarry Smith 
9280971522cSBarry Smith   PetscFunctionBegin;
9290700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
930db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
931db4c96c1SJed Brown   if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
932db4c96c1SJed Brown   PetscValidIntPointer(fields,3);
9330971522cSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr);
9340971522cSBarry Smith   if (f) {
935db4c96c1SJed Brown     ierr = (*f)(pc,splitname,n,fields);CHKERRQ(ierr);
9360971522cSBarry Smith   }
9370971522cSBarry Smith   PetscFunctionReturn(0);
9380971522cSBarry Smith }
9390971522cSBarry Smith 
9400971522cSBarry Smith #undef __FUNCT__
941704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
942704ba839SBarry Smith /*@
943704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
944704ba839SBarry Smith 
945*ad4df100SBarry Smith     Logically Collective on PC
946704ba839SBarry Smith 
947704ba839SBarry Smith     Input Parameters:
948704ba839SBarry Smith +   pc  - the preconditioner context
949db4c96c1SJed Brown .   splitname - name of this split
950db4c96c1SJed Brown -   is - the index set that defines the vector elements in this field
951704ba839SBarry Smith 
952d32f9abdSBarry Smith 
953a6ffb8dbSJed Brown     Notes:
954a6ffb8dbSJed Brown     Use PCFieldSplitSetFields(), for fields defined by strided types.
955a6ffb8dbSJed Brown 
956db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
957db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
958d32f9abdSBarry Smith 
959704ba839SBarry Smith     Level: intermediate
960704ba839SBarry Smith 
961704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
962704ba839SBarry Smith 
963704ba839SBarry Smith @*/
964db4c96c1SJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
965704ba839SBarry Smith {
966db4c96c1SJed Brown   PetscErrorCode ierr,(*f)(PC,const char[],IS);
967704ba839SBarry Smith 
968704ba839SBarry Smith   PetscFunctionBegin;
9690700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
970db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
971db4c96c1SJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,3);
972704ba839SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr);
973704ba839SBarry Smith   if (f) {
974db4c96c1SJed Brown     ierr = (*f)(pc,splitname,is);CHKERRQ(ierr);
975704ba839SBarry Smith   }
976704ba839SBarry Smith   PetscFunctionReturn(0);
977704ba839SBarry Smith }
978704ba839SBarry Smith 
979704ba839SBarry Smith #undef __FUNCT__
98051f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
98151f519a2SBarry Smith /*@
98251f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
98351f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
98451f519a2SBarry Smith 
985*ad4df100SBarry Smith     Logically Collective on PC
98651f519a2SBarry Smith 
98751f519a2SBarry Smith     Input Parameters:
98851f519a2SBarry Smith +   pc  - the preconditioner context
98951f519a2SBarry Smith -   bs - the block size
99051f519a2SBarry Smith 
99151f519a2SBarry Smith     Level: intermediate
99251f519a2SBarry Smith 
99351f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
99451f519a2SBarry Smith 
99551f519a2SBarry Smith @*/
99651f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
99751f519a2SBarry Smith {
99851f519a2SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt);
99951f519a2SBarry Smith 
100051f519a2SBarry Smith   PetscFunctionBegin;
10010700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
100251f519a2SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr);
100351f519a2SBarry Smith   if (f) {
100451f519a2SBarry Smith     ierr = (*f)(pc,bs);CHKERRQ(ierr);
100551f519a2SBarry Smith   }
100651f519a2SBarry Smith   PetscFunctionReturn(0);
100751f519a2SBarry Smith }
100851f519a2SBarry Smith 
100951f519a2SBarry Smith #undef __FUNCT__
101069a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
10110971522cSBarry Smith /*@C
101269a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
10130971522cSBarry Smith 
101469a612a9SBarry Smith    Collective on KSP
10150971522cSBarry Smith 
10160971522cSBarry Smith    Input Parameter:
10170971522cSBarry Smith .  pc - the preconditioner context
10180971522cSBarry Smith 
10190971522cSBarry Smith    Output Parameters:
10200971522cSBarry Smith +  n - the number of split
102169a612a9SBarry Smith -  pc - the array of KSP contexts
10220971522cSBarry Smith 
10230971522cSBarry Smith    Note:
1024d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
1025d32f9abdSBarry Smith    (not the KSP just the array that contains them).
10260971522cSBarry Smith 
102769a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
10280971522cSBarry Smith 
10290971522cSBarry Smith    Level: advanced
10300971522cSBarry Smith 
10310971522cSBarry Smith .seealso: PCFIELDSPLIT
10320971522cSBarry Smith @*/
1033dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
10340971522cSBarry Smith {
103569a612a9SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **);
10360971522cSBarry Smith 
10370971522cSBarry Smith   PetscFunctionBegin;
10380700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10390971522cSBarry Smith   PetscValidIntPointer(n,2);
104069a612a9SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr);
10410971522cSBarry Smith   if (f) {
104269a612a9SBarry Smith     ierr = (*f)(pc,n,subksp);CHKERRQ(ierr);
1043e7e72b3dSBarry Smith   } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC");
10440971522cSBarry Smith   PetscFunctionReturn(0);
10450971522cSBarry Smith }
10460971522cSBarry Smith 
1047e69d4d44SBarry Smith #undef __FUNCT__
1048e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition"
1049e69d4d44SBarry Smith /*@
1050e69d4d44SBarry Smith     PCFieldSplitSchurPrecondition -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
1051e69d4d44SBarry Smith       D matrix. Otherwise no preconditioner is used.
1052e69d4d44SBarry Smith 
1053e69d4d44SBarry Smith     Collective on PC
1054e69d4d44SBarry Smith 
1055e69d4d44SBarry Smith     Input Parameters:
1056e69d4d44SBarry Smith +   pc  - the preconditioner context
1057084e4875SJed Brown .   ptype - which matrix to use for preconditioning the Schur complement
1058084e4875SJed Brown -   userpre - matrix to use for preconditioning, or PETSC_NULL
1059084e4875SJed Brown 
1060084e4875SJed Brown     Notes:
1061084e4875SJed Brown     The default is to use the block on the diagonal of the preconditioning matrix.  This is D, in the (1,1) position.
1062084e4875SJed Brown     There are currently no preconditioners that work directly with the Schur complement so setting
1063084e4875SJed Brown     PC_FIELDSPLIT_SCHUR_PRE_SELF is observationally equivalent to -fieldsplit_1_pc_type none.
1064e69d4d44SBarry Smith 
1065e69d4d44SBarry Smith     Options Database:
1066084e4875SJed Brown .     -pc_fieldsplit_schur_precondition <self,user,diag> default is diag
1067e69d4d44SBarry Smith 
1068e69d4d44SBarry Smith     Level: intermediate
1069e69d4d44SBarry Smith 
1070084e4875SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType
1071e69d4d44SBarry Smith 
1072e69d4d44SBarry Smith @*/
1073084e4875SJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1074e69d4d44SBarry Smith {
1075084e4875SJed Brown   PetscErrorCode ierr,(*f)(PC,PCFieldSplitSchurPreType,Mat);
1076e69d4d44SBarry Smith 
1077e69d4d44SBarry Smith   PetscFunctionBegin;
10780700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1079e69d4d44SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",(void (**)(void))&f);CHKERRQ(ierr);
1080e69d4d44SBarry Smith   if (f) {
1081084e4875SJed Brown     ierr = (*f)(pc,ptype,pre);CHKERRQ(ierr);
1082e69d4d44SBarry Smith   }
1083e69d4d44SBarry Smith   PetscFunctionReturn(0);
1084e69d4d44SBarry Smith }
1085e69d4d44SBarry Smith 
1086e69d4d44SBarry Smith EXTERN_C_BEGIN
1087e69d4d44SBarry Smith #undef __FUNCT__
1088e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
1089084e4875SJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1090e69d4d44SBarry Smith {
1091e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1092084e4875SJed Brown   PetscErrorCode  ierr;
1093e69d4d44SBarry Smith 
1094e69d4d44SBarry Smith   PetscFunctionBegin;
1095084e4875SJed Brown   jac->schurpre = ptype;
1096084e4875SJed Brown   if (pre) {
1097084e4875SJed Brown     if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);}
1098084e4875SJed Brown     jac->schur_user = pre;
1099084e4875SJed Brown     ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1100084e4875SJed Brown   }
1101e69d4d44SBarry Smith   PetscFunctionReturn(0);
1102e69d4d44SBarry Smith }
1103e69d4d44SBarry Smith EXTERN_C_END
1104e69d4d44SBarry Smith 
110530ad9308SMatthew Knepley #undef __FUNCT__
110630ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
110730ad9308SMatthew Knepley /*@C
110830ad9308SMatthew Knepley    PCFieldSplitGetSchurBlocks - Gets the all matrix blocks for the Schur complement
110930ad9308SMatthew Knepley 
111030ad9308SMatthew Knepley    Collective on KSP
111130ad9308SMatthew Knepley 
111230ad9308SMatthew Knepley    Input Parameter:
111330ad9308SMatthew Knepley .  pc - the preconditioner context
111430ad9308SMatthew Knepley 
111530ad9308SMatthew Knepley    Output Parameters:
111630ad9308SMatthew Knepley +  A - the (0,0) block
111730ad9308SMatthew Knepley .  B - the (0,1) block
111830ad9308SMatthew Knepley .  C - the (1,0) block
111930ad9308SMatthew Knepley -  D - the (1,1) block
112030ad9308SMatthew Knepley 
112130ad9308SMatthew Knepley    Level: advanced
112230ad9308SMatthew Knepley 
112330ad9308SMatthew Knepley .seealso: PCFIELDSPLIT
112430ad9308SMatthew Knepley @*/
112530ad9308SMatthew Knepley PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSchurBlocks(PC pc,Mat *A,Mat *B,Mat *C, Mat *D)
112630ad9308SMatthew Knepley {
112730ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit *) pc->data;
112830ad9308SMatthew Knepley 
112930ad9308SMatthew Knepley   PetscFunctionBegin;
11300700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1131c1235816SBarry Smith   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
113230ad9308SMatthew Knepley   if (A) *A = jac->pmat[0];
113330ad9308SMatthew Knepley   if (B) *B = jac->B;
113430ad9308SMatthew Knepley   if (C) *C = jac->C;
113530ad9308SMatthew Knepley   if (D) *D = jac->pmat[1];
113630ad9308SMatthew Knepley   PetscFunctionReturn(0);
113730ad9308SMatthew Knepley }
113830ad9308SMatthew Knepley 
113979416396SBarry Smith EXTERN_C_BEGIN
114079416396SBarry Smith #undef __FUNCT__
114179416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
1142dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
114379416396SBarry Smith {
114479416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1145e69d4d44SBarry Smith   PetscErrorCode ierr;
114679416396SBarry Smith 
114779416396SBarry Smith   PetscFunctionBegin;
114879416396SBarry Smith   jac->type = type;
11493b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
11503b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
11513b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
1152e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1153e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1154e69d4d44SBarry Smith 
11553b224e63SBarry Smith   } else {
11563b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
11573b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
1158e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
11599e7d6b0aSBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr);
11603b224e63SBarry Smith   }
116179416396SBarry Smith   PetscFunctionReturn(0);
116279416396SBarry Smith }
116379416396SBarry Smith EXTERN_C_END
116479416396SBarry Smith 
116551f519a2SBarry Smith EXTERN_C_BEGIN
116651f519a2SBarry Smith #undef __FUNCT__
116751f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
116851f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
116951f519a2SBarry Smith {
117051f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
117151f519a2SBarry Smith 
117251f519a2SBarry Smith   PetscFunctionBegin;
117365e19b50SBarry Smith   if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
117465e19b50SBarry 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);
117551f519a2SBarry Smith   jac->bs = bs;
117651f519a2SBarry Smith   PetscFunctionReturn(0);
117751f519a2SBarry Smith }
117851f519a2SBarry Smith EXTERN_C_END
117951f519a2SBarry Smith 
118079416396SBarry Smith #undef __FUNCT__
118179416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
1182bc08b0f1SBarry Smith /*@
118379416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
118479416396SBarry Smith 
118579416396SBarry Smith    Collective on PC
118679416396SBarry Smith 
118779416396SBarry Smith    Input Parameter:
118879416396SBarry Smith .  pc - the preconditioner context
118981540f2fSBarry Smith .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
119079416396SBarry Smith 
119179416396SBarry Smith    Options Database Key:
1192a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
119379416396SBarry Smith 
119479416396SBarry Smith    Level: Developer
119579416396SBarry Smith 
119679416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
119779416396SBarry Smith 
119879416396SBarry Smith .seealso: PCCompositeSetType()
119979416396SBarry Smith 
120079416396SBarry Smith @*/
1201dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type)
120279416396SBarry Smith {
120379416396SBarry Smith   PetscErrorCode ierr,(*f)(PC,PCCompositeType);
120479416396SBarry Smith 
120579416396SBarry Smith   PetscFunctionBegin;
12060700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
120779416396SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
120879416396SBarry Smith   if (f) {
120979416396SBarry Smith     ierr = (*f)(pc,type);CHKERRQ(ierr);
121079416396SBarry Smith   }
121179416396SBarry Smith   PetscFunctionReturn(0);
121279416396SBarry Smith }
121379416396SBarry Smith 
12140971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
12150971522cSBarry Smith /*MC
1216a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
12170971522cSBarry Smith                   fields or groups of fields
12180971522cSBarry Smith 
12190971522cSBarry Smith 
1220edf189efSBarry Smith      To set options on the solvers for each block append -fieldsplit_ to all the PC
1221edf189efSBarry Smith         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
12220971522cSBarry Smith 
1223a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
122469a612a9SBarry Smith          and set the options directly on the resulting KSP object
12250971522cSBarry Smith 
12260971522cSBarry Smith    Level: intermediate
12270971522cSBarry Smith 
122879416396SBarry Smith    Options Database Keys:
122981540f2fSBarry Smith +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
123081540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
123181540f2fSBarry Smith                               been supplied explicitly by -pc_fieldsplit_%d_fields
123281540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
123381540f2fSBarry Smith .    -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative>
1234e69d4d44SBarry Smith .    -pc_fieldsplit_schur_precondition <true,false> default is true
123579416396SBarry Smith 
1236edf189efSBarry Smith -    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
12373b224e63SBarry Smith      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
12383b224e63SBarry Smith 
12393b224e63SBarry Smith 
1240d32f9abdSBarry Smith    Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1241d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
1242d32f9abdSBarry Smith 
1243d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
1244d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1245d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
1246d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
1247d32f9abdSBarry Smith 
1248d32f9abdSBarry Smith       Currently for the multiplicative version, the updated residual needed for the next field
1249d32f9abdSBarry Smith      solve is computed via a matrix vector product over the entire array. An optimization would be
1250d32f9abdSBarry Smith      to update the residual only for the part of the right hand side associated with the next field
1251d32f9abdSBarry Smith      solve. (This would involve more MatGetSubMatrix() calls or some other mechanism to compute the
1252d32f9abdSBarry Smith      part of the matrix needed to just update part of the residual).
1253d32f9abdSBarry Smith 
1254e69d4d44SBarry Smith      For the Schur complement preconditioner if J = ( A B )
1255e69d4d44SBarry Smith                                                     ( C D )
1256e69d4d44SBarry Smith      the preconditioner is
1257e69d4d44SBarry Smith               (I   -B inv(A)) ( inv(A)   0    ) (I         0  )
1258e69d4d44SBarry Smith               (0    I       ) (   0    inv(S) ) (-C inv(A) I  )
1259edf189efSBarry Smith      where the action of inv(A) is applied using the KSP solver with prefix -fieldsplit_0_. The action of
1260e69d4d44SBarry Smith      inv(S) is computed using the KSP solver with prefix -schur_. For PCFieldSplitGetKSP() when field number is
1261edf189efSBarry Smith      0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
1262e69d4d44SBarry Smith      D is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
1263e69d4d44SBarry Smith      option.
1264e69d4d44SBarry Smith 
1265edf189efSBarry Smith      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
1266edf189efSBarry Smith      is used automatically for a second block.
1267edf189efSBarry Smith 
1268a541d17aSBarry Smith    Concepts: physics based preconditioners, block preconditioners
12690971522cSBarry Smith 
1270a541d17aSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners
1271e69d4d44SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
12720971522cSBarry Smith M*/
12730971522cSBarry Smith 
12740971522cSBarry Smith EXTERN_C_BEGIN
12750971522cSBarry Smith #undef __FUNCT__
12760971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
1277dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc)
12780971522cSBarry Smith {
12790971522cSBarry Smith   PetscErrorCode ierr;
12800971522cSBarry Smith   PC_FieldSplit  *jac;
12810971522cSBarry Smith 
12820971522cSBarry Smith   PetscFunctionBegin;
128338f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
12840971522cSBarry Smith   jac->bs        = -1;
12850971522cSBarry Smith   jac->nsplits   = 0;
12863e197d65SBarry Smith   jac->type      = PC_COMPOSITE_MULTIPLICATIVE;
1287e6cab6aaSJed Brown   jac->schurpre  = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
1288c5d2311dSJed Brown   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL;
128951f519a2SBarry Smith 
12900971522cSBarry Smith   pc->data     = (void*)jac;
12910971522cSBarry Smith 
12920971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
1293421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
12940971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
12950971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
12960971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
12970971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
12980971522cSBarry Smith   pc->ops->applyrichardson   = 0;
12990971522cSBarry Smith 
130069a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
130169a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
13020971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
13030971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1304704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
1305704ba839SBarry Smith                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
130679416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
130779416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
130851f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
130951f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
13100971522cSBarry Smith   PetscFunctionReturn(0);
13110971522cSBarry Smith }
13120971522cSBarry Smith EXTERN_C_END
13130971522cSBarry Smith 
13140971522cSBarry Smith 
1315a541d17aSBarry Smith /*MC
1316a541d17aSBarry Smith    Block_Preconditioners - PETSc provides support for a variety of "block preconditioners", this provides an
1317a541d17aSBarry Smith           overview of these methods.
1318a541d17aSBarry Smith 
1319a541d17aSBarry Smith       Consider the solution to ( A B ) (x_1)  =  (b_1)
1320a541d17aSBarry Smith                                ( C D ) (x_2)     (b_2)
1321a541d17aSBarry Smith 
1322a541d17aSBarry Smith       Important special cases, the Stokes equation: C = B' and D = 0  (A   B) (x_1) = (b_1)
1323a541d17aSBarry Smith                                                                        B'  0) (x_2)   (b_2)
1324a541d17aSBarry Smith 
1325a541d17aSBarry Smith       One of the goals of the PCFieldSplit preconditioner in PETSc is to provide a variety of preconditioners
1326a541d17aSBarry Smith       for this block system.
1327a541d17aSBarry Smith 
1328a541d17aSBarry Smith       Consider an additional matrix (Ap  Bp)
1329a541d17aSBarry Smith                                     (Cp  Dp) where some or all of the entries may be the same as
1330a541d17aSBarry Smith       in the original matrix (for example Ap == A).
1331a541d17aSBarry Smith 
1332a541d17aSBarry Smith       In the following, A^ denotes the approximate application of the inverse of A, possibly using Ap in the
1333a541d17aSBarry Smith       approximation. In PETSc this simply means one has called KSPSetOperators(ksp,A,Ap,...) or KSPSetOperators(ksp,Ap,Ap,...)
1334a541d17aSBarry Smith 
1335a541d17aSBarry Smith       Block Jacobi:   x_1 = A^ b_1
1336a541d17aSBarry Smith                       x_2 = D^ b_2
1337a541d17aSBarry Smith 
1338a541d17aSBarry Smith       Lower block Gauss-Seidel:   x_1 = A^ b_1
1339a541d17aSBarry Smith                             x_2 = D^ (b_2 - C x_1)       variant x_2 = D^ (b_2 - Cp x_1)
1340a541d17aSBarry Smith 
1341a541d17aSBarry Smith       Symmetric Gauss-Seidel:  x_1 = x_1 + A^(b_1 - A x_1 - B x_2)    variant  x_1 = x_1 + A^(b_1 - Ap x_1 - Bp x_2)
1342a541d17aSBarry Smith           Interestingly this form is not actually a symmetric matrix, the symmetric version is
1343a541d17aSBarry Smith                               x_1 = A^(b_1 - B x_2)      variant x_1 = A^(b_1 - Bp x_2)
1344a541d17aSBarry Smith 
13450bc0a719SBarry Smith    Level: intermediate
13460bc0a719SBarry Smith 
1347a541d17aSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCFIELDSPLIT
1348a541d17aSBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
1349a541d17aSBarry Smith M*/
1350