xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision d0af7cd330563f97edd88b8a891ae212dc5b714a)
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;
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;
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;
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;
76ace3abfcSBarry Smith   PetscBool         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;
117ace3abfcSBarry Smith   PetscBool         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     }
145435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 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);
154435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \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;
177ace3abfcSBarry Smith   PetscBool      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;
2086ce1633cSBarry Smith   PetscBool         flg = PETSC_FALSE,stokes = PETSC_FALSE;
2096c924f48SJed Brown   PetscInt          i;
2100971522cSBarry Smith 
2110971522cSBarry Smith   PetscFunctionBegin;
212d32f9abdSBarry Smith   if (!ilink) {
213*d0af7cd3SBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
214*d0af7cd3SBarry Smith     if (pc->dm && !stokes) {
2158b8307b2SJed Brown       PetscBool dmcomposite;
2168b8307b2SJed Brown       ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr);
2178b8307b2SJed Brown       if (dmcomposite) {
2188b8307b2SJed Brown         PetscInt nDM;
2198b8307b2SJed Brown         IS       *fields;
2208b8307b2SJed Brown         ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
2218b8307b2SJed Brown         ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr);
2228b8307b2SJed Brown         ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr);
2238b8307b2SJed Brown         for (i=0; i<nDM; i++) {
2248b8307b2SJed Brown           char splitname[8];
2258b8307b2SJed Brown           ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
2268b8307b2SJed Brown           ierr = PCFieldSplitSetIS(pc,splitname,fields[i]);CHKERRQ(ierr);
227fcfd50ebSBarry Smith           ierr = ISDestroy(&fields[i]);CHKERRQ(ierr);
2288b8307b2SJed Brown         }
2298b8307b2SJed Brown         ierr = PetscFree(fields);CHKERRQ(ierr);
2308b8307b2SJed Brown       }
23166ffff09SJed Brown     } else {
232521d7252SBarry Smith       if (jac->bs <= 0) {
233704ba839SBarry Smith         if (pc->pmat) {
234521d7252SBarry Smith           ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
235704ba839SBarry Smith         } else {
236704ba839SBarry Smith           jac->bs = 1;
237704ba839SBarry Smith         }
238521d7252SBarry Smith       }
239d32f9abdSBarry Smith 
240acfcf0e5SJed Brown       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
2416ce1633cSBarry Smith       if (stokes) {
2426ce1633cSBarry Smith         IS       zerodiags,rest;
2436ce1633cSBarry Smith         PetscInt nmin,nmax;
2446ce1633cSBarry Smith 
2456ce1633cSBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
2466ce1633cSBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
2476ce1633cSBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
2486ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
2496ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr);
250fcfd50ebSBarry Smith         ierr = ISDestroy(&zerodiags);CHKERRQ(ierr);
251fcfd50ebSBarry Smith         ierr = ISDestroy(&rest);CHKERRQ(ierr);
2526ce1633cSBarry Smith       } else {
253d32f9abdSBarry Smith         if (!flg) {
254d32f9abdSBarry Smith           /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
255d32f9abdSBarry Smith            then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
2566c924f48SJed Brown           ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
2576c924f48SJed Brown           if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
258d32f9abdSBarry Smith         }
2596c924f48SJed Brown         if (flg || !jac->splitdefined) {
260d32f9abdSBarry Smith           ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
261db4c96c1SJed Brown           for (i=0; i<jac->bs; i++) {
2626c924f48SJed Brown             char splitname[8];
2636c924f48SJed Brown             ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
264db4c96c1SJed Brown             ierr = PCFieldSplitSetFields(pc,splitname,1,&i);CHKERRQ(ierr);
26579416396SBarry Smith           }
26697bbdb24SBarry Smith           jac->defaultsplit = PETSC_TRUE;
267521d7252SBarry Smith         }
26866ffff09SJed Brown       }
2696ce1633cSBarry Smith     }
270edf189efSBarry Smith   } else if (jac->nsplits == 1) {
271edf189efSBarry Smith     if (ilink->is) {
272edf189efSBarry Smith       IS       is2;
273edf189efSBarry Smith       PetscInt nmin,nmax;
274edf189efSBarry Smith 
275edf189efSBarry Smith       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
276edf189efSBarry Smith       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
277db4c96c1SJed Brown       ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr);
278fcfd50ebSBarry Smith       ierr = ISDestroy(&is2);CHKERRQ(ierr);
279e7e72b3dSBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
280*d0af7cd3SBarry Smith   } else if (!pc->setupcalled) {
281*d0af7cd3SBarry Smith     /* PCReset() has been called on this PC, ilink exists but all data structures in it must be rebuilt
282*d0af7cd3SBarry Smith        This is basically the !ilink portion of code above copied from above and the allocation of the ilinks removed
283*d0af7cd3SBarry Smith        since they already exist. This should be totally rewritten */
284*d0af7cd3SBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
285*d0af7cd3SBarry Smith     if (pc->dm && !stokes) {
286*d0af7cd3SBarry Smith       PetscBool dmcomposite;
287*d0af7cd3SBarry Smith       ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr);
288*d0af7cd3SBarry Smith       if (dmcomposite) {
289*d0af7cd3SBarry Smith         PetscInt nDM;
290*d0af7cd3SBarry Smith         IS       *fields;
291*d0af7cd3SBarry Smith         ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
292*d0af7cd3SBarry Smith         ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr);
293*d0af7cd3SBarry Smith         ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr);
294*d0af7cd3SBarry Smith         for (i=0; i<nDM; i++) {
295*d0af7cd3SBarry Smith           ilink->is = fields[i];
296*d0af7cd3SBarry Smith           ilink     = ilink->next;
297edf189efSBarry Smith         }
298*d0af7cd3SBarry Smith         ierr = PetscFree(fields);CHKERRQ(ierr);
299*d0af7cd3SBarry Smith       }
300*d0af7cd3SBarry Smith     } else {
301*d0af7cd3SBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
302*d0af7cd3SBarry Smith       if (stokes) {
303*d0af7cd3SBarry Smith         IS       zerodiags,rest;
304*d0af7cd3SBarry Smith         PetscInt nmin,nmax;
305*d0af7cd3SBarry Smith 
306*d0af7cd3SBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
307*d0af7cd3SBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
308*d0af7cd3SBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
309*d0af7cd3SBarry Smith         ilink->is       = rest;
310*d0af7cd3SBarry Smith         ilink->next->is = zerodiags;
311*d0af7cd3SBarry Smith       } else {
312*d0af7cd3SBarry Smith         SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used");
313*d0af7cd3SBarry Smith       }
314*d0af7cd3SBarry Smith     }
315*d0af7cd3SBarry Smith   }
316*d0af7cd3SBarry Smith 
317e7e72b3dSBarry Smith   if (jac->nsplits < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields");
31869a612a9SBarry Smith   PetscFunctionReturn(0);
31969a612a9SBarry Smith }
32069a612a9SBarry Smith 
32169a612a9SBarry Smith #undef __FUNCT__
32269a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
32369a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
32469a612a9SBarry Smith {
32569a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
32669a612a9SBarry Smith   PetscErrorCode    ierr;
3275a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
328cf502942SBarry Smith   PetscInt          i,nsplit,ccsize;
32969a612a9SBarry Smith   MatStructure      flag = pc->flag;
330ace3abfcSBarry Smith   PetscBool         sorted;
33169a612a9SBarry Smith 
33269a612a9SBarry Smith   PetscFunctionBegin;
33369a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
33497bbdb24SBarry Smith   nsplit = jac->nsplits;
3355a9f2f41SSatish Balay   ilink  = jac->head;
33697bbdb24SBarry Smith 
33797bbdb24SBarry Smith   /* get the matrices for each split */
338704ba839SBarry Smith   if (!jac->issetup) {
3391b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
34097bbdb24SBarry Smith 
341704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
342704ba839SBarry Smith 
343704ba839SBarry Smith     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
34451f519a2SBarry Smith     bs     = jac->bs;
34597bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
346cf502942SBarry Smith     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
3471b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
3481b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
3491b9fc7fcSBarry Smith       if (jac->defaultsplit) {
350704ba839SBarry Smith         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
351704ba839SBarry Smith       } else if (!ilink->is) {
352ccb205f8SBarry Smith         if (ilink->nfields > 1) {
3535a9f2f41SSatish Balay           PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
3545a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
3551b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
3561b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
3571b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
35897bbdb24SBarry Smith             }
35997bbdb24SBarry Smith           }
360d67e408aSBarry Smith           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr);
361ccb205f8SBarry Smith         } else {
362704ba839SBarry Smith           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
363ccb205f8SBarry Smith         }
3643e197d65SBarry Smith       }
365704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
366e32f2f54SBarry Smith       if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
367704ba839SBarry Smith       ilink = ilink->next;
3681b9fc7fcSBarry Smith     }
3691b9fc7fcSBarry Smith   }
3701b9fc7fcSBarry Smith 
371704ba839SBarry Smith   ilink  = jac->head;
37297bbdb24SBarry Smith   if (!jac->pmat) {
373cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
374cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
3754aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
376704ba839SBarry Smith       ilink = ilink->next;
377cf502942SBarry Smith     }
37897bbdb24SBarry Smith   } else {
379cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
3804aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
381704ba839SBarry Smith       ilink = ilink->next;
382cf502942SBarry Smith     }
38397bbdb24SBarry Smith   }
384519d70e2SJed Brown   if (jac->realdiagonal) {
385519d70e2SJed Brown     ilink = jac->head;
386519d70e2SJed Brown     if (!jac->mat) {
387519d70e2SJed Brown       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr);
388519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
389519d70e2SJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
390519d70e2SJed Brown         ilink = ilink->next;
391519d70e2SJed Brown       }
392519d70e2SJed Brown     } else {
393519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
394966e74d7SJed Brown         if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);}
395519d70e2SJed Brown         ilink = ilink->next;
396519d70e2SJed Brown       }
397519d70e2SJed Brown     }
398519d70e2SJed Brown   } else {
399519d70e2SJed Brown     jac->mat = jac->pmat;
400519d70e2SJed Brown   }
40197bbdb24SBarry Smith 
4026c8605c2SJed Brown   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
40368dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
40468dd23aaSBarry Smith     ilink  = jac->head;
40568dd23aaSBarry Smith     if (!jac->Afield) {
40668dd23aaSBarry Smith       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
40768dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
4084aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
40968dd23aaSBarry Smith         ilink = ilink->next;
41068dd23aaSBarry Smith       }
41168dd23aaSBarry Smith     } else {
41268dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
4134aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
41468dd23aaSBarry Smith         ilink = ilink->next;
41568dd23aaSBarry Smith       }
41668dd23aaSBarry Smith     }
41768dd23aaSBarry Smith   }
41868dd23aaSBarry Smith 
4193b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
4203b224e63SBarry Smith     IS       ccis;
4214aa3045dSJed Brown     PetscInt rstart,rend;
422e7e72b3dSBarry Smith     if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
42368dd23aaSBarry Smith 
424e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
425e6cab6aaSJed Brown     ierr  = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
426e6cab6aaSJed Brown 
4273b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
4283b224e63SBarry Smith     if (jac->schur) {
4293b224e63SBarry Smith       ilink = jac->head;
4304aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4314aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
432fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4333b224e63SBarry Smith       ilink = ilink->next;
4344aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4354aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
436fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
437519d70e2SJed Brown       ierr  = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr);
438084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr);
4393b224e63SBarry Smith 
4403b224e63SBarry Smith      } else {
4411cee3971SBarry Smith       KSP ksp;
4426c924f48SJed Brown       char schurprefix[256];
4433b224e63SBarry Smith 
444a04f6461SBarry Smith       /* extract the A01 and A10 matrices */
4453b224e63SBarry Smith       ilink = jac->head;
4464aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4474aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
448fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4493b224e63SBarry Smith       ilink = ilink->next;
4504aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4514aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
452fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4537d50072eSJed Brown       /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */
454519d70e2SJed Brown       ierr  = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr);
455a04f6461SBarry Smith       /* set tabbing and options prefix of KSP inside the MatSchur */
4561cee3971SBarry Smith       ierr  = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
457e69d4d44SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr);
458a04f6461SBarry Smith       ierr  = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",jac->head->splitname);CHKERRQ(ierr);
459a04f6461SBarry Smith       ierr  = KSPSetOptionsPrefix(ksp,schurprefix);CHKERRQ(ierr);
4603b224e63SBarry Smith       ierr  = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
4613b224e63SBarry Smith 
4623b224e63SBarry Smith       ierr  = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr);
4639005cf84SBarry Smith       ierr  = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr);
4641cee3971SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
465084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
466084e4875SJed Brown       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
467084e4875SJed Brown         PC pc;
468cd5f4a64SJed Brown         ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr);
469084e4875SJed Brown         ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
470084e4875SJed Brown         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
471e69d4d44SBarry Smith       }
4726c924f48SJed Brown       ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
4736c924f48SJed Brown       ierr  = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr);
4743b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
4753b224e63SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
4763b224e63SBarry Smith 
4773b224e63SBarry Smith       ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr);
4783b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr);
4793b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr);
4803b224e63SBarry Smith       ilink = jac->head;
4813b224e63SBarry Smith       ilink->x = jac->x[0]; ilink->y = jac->y[0];
4823b224e63SBarry Smith       ilink = ilink->next;
4833b224e63SBarry Smith       ilink->x = jac->x[1]; ilink->y = jac->y[1];
4843b224e63SBarry Smith     }
4853b224e63SBarry Smith   } else {
48697bbdb24SBarry Smith     /* set up the individual PCs */
48797bbdb24SBarry Smith     i    = 0;
4885a9f2f41SSatish Balay     ilink = jac->head;
4895a9f2f41SSatish Balay     while (ilink) {
490519d70e2SJed Brown       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr);
4913b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
4925a9f2f41SSatish Balay       ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);
4935a9f2f41SSatish Balay       ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
49497bbdb24SBarry Smith       i++;
4955a9f2f41SSatish Balay       ilink = ilink->next;
4960971522cSBarry Smith     }
49797bbdb24SBarry Smith 
49897bbdb24SBarry Smith     /* create work vectors for each split */
4991b9fc7fcSBarry Smith     if (!jac->x) {
50097bbdb24SBarry Smith       ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
5015a9f2f41SSatish Balay       ilink = jac->head;
50297bbdb24SBarry Smith       for (i=0; i<nsplit; i++) {
503906ed7ccSBarry Smith         Vec *vl,*vr;
504906ed7ccSBarry Smith 
505906ed7ccSBarry Smith         ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
506906ed7ccSBarry Smith         ilink->x  = *vr;
507906ed7ccSBarry Smith         ilink->y  = *vl;
508906ed7ccSBarry Smith         ierr      = PetscFree(vr);CHKERRQ(ierr);
509906ed7ccSBarry Smith         ierr      = PetscFree(vl);CHKERRQ(ierr);
5105a9f2f41SSatish Balay         jac->x[i] = ilink->x;
5115a9f2f41SSatish Balay         jac->y[i] = ilink->y;
5125a9f2f41SSatish Balay         ilink     = ilink->next;
51397bbdb24SBarry Smith       }
5143b224e63SBarry Smith     }
5153b224e63SBarry Smith   }
5163b224e63SBarry Smith 
5173b224e63SBarry Smith 
518d0f46423SBarry Smith   if (!jac->head->sctx) {
5193b224e63SBarry Smith     Vec xtmp;
5203b224e63SBarry Smith 
52179416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
5221b9fc7fcSBarry Smith 
5235a9f2f41SSatish Balay     ilink = jac->head;
5241b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
5251b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
526704ba839SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
5275a9f2f41SSatish Balay       ilink = ilink->next;
5281b9fc7fcSBarry Smith     }
5296bf464f9SBarry Smith     ierr = VecDestroy(&xtmp);CHKERRQ(ierr);
5301b9fc7fcSBarry Smith   }
5310971522cSBarry Smith   PetscFunctionReturn(0);
5320971522cSBarry Smith }
5330971522cSBarry Smith 
5345a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
535ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
536ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
5375a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
538ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
539ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
54079416396SBarry Smith 
5410971522cSBarry Smith #undef __FUNCT__
5423b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur"
5433b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
5443b224e63SBarry Smith {
5453b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
5463b224e63SBarry Smith   PetscErrorCode    ierr;
5473b224e63SBarry Smith   KSP               ksp;
5483b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
5493b224e63SBarry Smith 
5503b224e63SBarry Smith   PetscFunctionBegin;
5513b224e63SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
5523b224e63SBarry Smith 
553c5d2311dSJed Brown   switch (jac->schurfactorization) {
554c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG:
555a04f6461SBarry Smith       /* [A00 0; 0 -S], positive definite, suitable for MINRES */
556c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
557c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
558c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
559c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
560c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
561c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
562c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
563c5d2311dSJed Brown       ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr);
564c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
565c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
566c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
567c5d2311dSJed Brown       break;
568c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER:
569a04f6461SBarry Smith       /* [A00 0; A10 S], suitable for left preconditioning */
570c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
571c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
572c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
573c5d2311dSJed Brown       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
574c5d2311dSJed Brown       ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr);
575c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);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,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
578c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
579c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
580c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
581c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
582c5d2311dSJed Brown       break;
583c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER:
584a04f6461SBarry Smith       /* [A00 A01; 0 S], suitable for right preconditioning */
585c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
586c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
587c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
588c5d2311dSJed Brown       ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
589c5d2311dSJed Brown       ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr);
590c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
591c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
592c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
593c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
594c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
595c5d2311dSJed Brown       ierr = VecScatterEnd(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       break;
598c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL:
599a04f6461SBarry 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 */
6003b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6013b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6023b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
6033b224e63SBarry Smith       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
6043b224e63SBarry Smith       ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
6053b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6063b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6073b224e63SBarry Smith 
6083b224e63SBarry Smith       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
6093b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6103b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6113b224e63SBarry Smith 
6123b224e63SBarry Smith       ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
6133b224e63SBarry Smith       ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
6143b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
6153b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6163b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
617c5d2311dSJed Brown   }
6183b224e63SBarry Smith   PetscFunctionReturn(0);
6193b224e63SBarry Smith }
6203b224e63SBarry Smith 
6213b224e63SBarry Smith #undef __FUNCT__
6220971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
6230971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
6240971522cSBarry Smith {
6250971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
6260971522cSBarry Smith   PetscErrorCode    ierr;
6275a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
628af93645bSJed Brown   PetscInt          cnt;
6290971522cSBarry Smith 
6300971522cSBarry Smith   PetscFunctionBegin;
63151f519a2SBarry Smith   CHKMEMQ;
63251f519a2SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
63351f519a2SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
63451f519a2SBarry Smith 
63579416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
6361b9fc7fcSBarry Smith     if (jac->defaultsplit) {
6370971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
6385a9f2f41SSatish Balay       while (ilink) {
6395a9f2f41SSatish Balay         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
6405a9f2f41SSatish Balay         ilink = ilink->next;
6410971522cSBarry Smith       }
6420971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
6431b9fc7fcSBarry Smith     } else {
644efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
6455a9f2f41SSatish Balay       while (ilink) {
6465a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
6475a9f2f41SSatish Balay         ilink = ilink->next;
6481b9fc7fcSBarry Smith       }
6491b9fc7fcSBarry Smith     }
65016913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
65179416396SBarry Smith     if (!jac->w1) {
65279416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
65379416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
65479416396SBarry Smith     }
655efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
6565a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
6573e197d65SBarry Smith     cnt = 1;
6585a9f2f41SSatish Balay     while (ilink->next) {
6595a9f2f41SSatish Balay       ilink = ilink->next;
6603e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
6613e197d65SBarry Smith       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
6623e197d65SBarry Smith       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
6633e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6643e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6653e197d65SBarry Smith       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
6663e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6673e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6683e197d65SBarry Smith     }
66951f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
67011755939SBarry Smith       cnt -= 2;
67151f519a2SBarry Smith       while (ilink->previous) {
67251f519a2SBarry Smith         ilink = ilink->previous;
67311755939SBarry Smith         /* compute the residual only over the part of the vector needed */
67411755939SBarry Smith         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
67511755939SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
67611755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
67711755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
67811755939SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
67911755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
68011755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
68151f519a2SBarry Smith       }
68211755939SBarry Smith     }
68365e19b50SBarry Smith   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
68451f519a2SBarry Smith   CHKMEMQ;
6850971522cSBarry Smith   PetscFunctionReturn(0);
6860971522cSBarry Smith }
6870971522cSBarry Smith 
688421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
689ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
690ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
691421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
692ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
693ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
694421e10b8SBarry Smith 
695421e10b8SBarry Smith #undef __FUNCT__
6968c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit"
697421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
698421e10b8SBarry Smith {
699421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
700421e10b8SBarry Smith   PetscErrorCode    ierr;
701421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
702421e10b8SBarry Smith 
703421e10b8SBarry Smith   PetscFunctionBegin;
704421e10b8SBarry Smith   CHKMEMQ;
705421e10b8SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
706421e10b8SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
707421e10b8SBarry Smith 
708421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
709421e10b8SBarry Smith     if (jac->defaultsplit) {
710421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
711421e10b8SBarry Smith       while (ilink) {
712421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
713421e10b8SBarry Smith 	ilink = ilink->next;
714421e10b8SBarry Smith       }
715421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
716421e10b8SBarry Smith     } else {
717421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
718421e10b8SBarry Smith       while (ilink) {
719421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
720421e10b8SBarry Smith 	ilink = ilink->next;
721421e10b8SBarry Smith       }
722421e10b8SBarry Smith     }
723421e10b8SBarry Smith   } else {
724421e10b8SBarry Smith     if (!jac->w1) {
725421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
726421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
727421e10b8SBarry Smith     }
728421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
729421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
730421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
731421e10b8SBarry Smith       while (ilink->next) {
732421e10b8SBarry Smith         ilink = ilink->next;
7339989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
734421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
735421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
736421e10b8SBarry Smith       }
737421e10b8SBarry Smith       while (ilink->previous) {
738421e10b8SBarry Smith         ilink = ilink->previous;
7399989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
740421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
741421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
742421e10b8SBarry Smith       }
743421e10b8SBarry Smith     } else {
744421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
745421e10b8SBarry Smith 	ilink = ilink->next;
746421e10b8SBarry Smith       }
747421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
748421e10b8SBarry Smith       while (ilink->previous) {
749421e10b8SBarry Smith 	ilink = ilink->previous;
7509989ab13SBarry Smith 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
751421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
752421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
753421e10b8SBarry Smith       }
754421e10b8SBarry Smith     }
755421e10b8SBarry Smith   }
756421e10b8SBarry Smith   CHKMEMQ;
757421e10b8SBarry Smith   PetscFunctionReturn(0);
758421e10b8SBarry Smith }
759421e10b8SBarry Smith 
7600971522cSBarry Smith #undef __FUNCT__
761574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit"
762574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc)
7630971522cSBarry Smith {
7640971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
7650971522cSBarry Smith   PetscErrorCode    ierr;
7665a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
7670971522cSBarry Smith 
7680971522cSBarry Smith   PetscFunctionBegin;
7695a9f2f41SSatish Balay   while (ilink) {
770574deadeSBarry Smith     ierr = KSPReset(ilink->ksp);CHKERRQ(ierr);
771fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->x);CHKERRQ(ierr);
772fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->y);CHKERRQ(ierr);
773fcfd50ebSBarry Smith     ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr);
774fcfd50ebSBarry Smith     ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
7755a9f2f41SSatish Balay     next = ilink->next;
7765a9f2f41SSatish Balay     ilink = next;
7770971522cSBarry Smith   }
77805b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
779574deadeSBarry Smith   if (jac->mat && jac->mat != jac->pmat) {
780574deadeSBarry Smith     ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);
781574deadeSBarry Smith   } else if (jac->mat) {
782574deadeSBarry Smith     jac->mat = PETSC_NULL;
783574deadeSBarry Smith   }
78497bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
78568dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
7866bf464f9SBarry Smith   ierr = VecDestroy(&jac->w1);CHKERRQ(ierr);
7876bf464f9SBarry Smith   ierr = VecDestroy(&jac->w2);CHKERRQ(ierr);
7886bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur);CHKERRQ(ierr);
7896bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
7906bf464f9SBarry Smith   ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr);
7916bf464f9SBarry Smith   ierr = MatDestroy(&jac->B);CHKERRQ(ierr);
7926bf464f9SBarry Smith   ierr = MatDestroy(&jac->C);CHKERRQ(ierr);
793574deadeSBarry Smith   PetscFunctionReturn(0);
794574deadeSBarry Smith }
795574deadeSBarry Smith 
796574deadeSBarry Smith #undef __FUNCT__
797574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
798574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
799574deadeSBarry Smith {
800574deadeSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
801574deadeSBarry Smith   PetscErrorCode    ierr;
802574deadeSBarry Smith   PC_FieldSplitLink ilink = jac->head,next;
803574deadeSBarry Smith 
804574deadeSBarry Smith   PetscFunctionBegin;
805574deadeSBarry Smith   ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr);
806574deadeSBarry Smith   while (ilink) {
8076bf464f9SBarry Smith     ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr);
808574deadeSBarry Smith     next = ilink->next;
809574deadeSBarry Smith     ierr = PetscFree(ilink->splitname);CHKERRQ(ierr);
810574deadeSBarry Smith     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
811574deadeSBarry Smith     ierr = PetscFree(ilink);CHKERRQ(ierr);
812574deadeSBarry Smith     ilink = next;
813574deadeSBarry Smith   }
814574deadeSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
815c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
8160971522cSBarry Smith   PetscFunctionReturn(0);
8170971522cSBarry Smith }
8180971522cSBarry Smith 
8190971522cSBarry Smith #undef __FUNCT__
8200971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
8210971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
8220971522cSBarry Smith {
8231b9fc7fcSBarry Smith   PetscErrorCode  ierr;
8246c924f48SJed Brown   PetscInt        bs;
825bc59fbc5SBarry Smith   PetscBool       flg,stokes = PETSC_FALSE;
8269dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
8273b224e63SBarry Smith   PCCompositeType ctype;
8281b9fc7fcSBarry Smith 
8290971522cSBarry Smith   PetscFunctionBegin;
8301b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
831acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr);
83251f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
83351f519a2SBarry Smith   if (flg) {
83451f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
83551f519a2SBarry Smith   }
836704ba839SBarry Smith 
837435f959eSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
838c0adfefeSBarry Smith   if (stokes) {
839c0adfefeSBarry Smith     ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
840c0adfefeSBarry Smith     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
841c0adfefeSBarry Smith   }
842c0adfefeSBarry Smith 
8433b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
8443b224e63SBarry Smith   if (flg) {
8453b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
8463b224e63SBarry Smith   }
847d32f9abdSBarry Smith 
848c30613efSMatthew Knepley   /* Only setup fields once */
849c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
850d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
851d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
8526c924f48SJed Brown     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
8536c924f48SJed Brown     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
854d32f9abdSBarry Smith   }
855c5d2311dSJed Brown   if (jac->type == PC_COMPOSITE_SCHUR) {
856c5d2311dSJed Brown     ierr = PetscOptionsEnum("-pc_fieldsplit_schur_factorization_type","Factorization to use","None",PCFieldSplitSchurFactorizationTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,PETSC_NULL);CHKERRQ(ierr);
857084e4875SJed 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);
858c5d2311dSJed Brown   }
8591b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
8600971522cSBarry Smith   PetscFunctionReturn(0);
8610971522cSBarry Smith }
8620971522cSBarry Smith 
8630971522cSBarry Smith /*------------------------------------------------------------------------------------*/
8640971522cSBarry Smith 
8650971522cSBarry Smith EXTERN_C_BEGIN
8660971522cSBarry Smith #undef __FUNCT__
8670971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
8687087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields)
8690971522cSBarry Smith {
87097bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
8710971522cSBarry Smith   PetscErrorCode    ierr;
8725a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
87369a612a9SBarry Smith   char              prefix[128];
87451f519a2SBarry Smith   PetscInt          i;
8750971522cSBarry Smith 
8760971522cSBarry Smith   PetscFunctionBegin;
8776c924f48SJed Brown   if (jac->splitdefined) {
8786c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
8796c924f48SJed Brown     PetscFunctionReturn(0);
8806c924f48SJed Brown   }
88151f519a2SBarry Smith   for (i=0; i<n; i++) {
882e32f2f54SBarry 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);
883e32f2f54SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
88451f519a2SBarry Smith   }
885704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
886a04f6461SBarry Smith   if (splitname) {
887db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
888a04f6461SBarry Smith   } else {
889a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
890a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
891a04f6461SBarry Smith   }
892704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
8935a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
8945a9f2f41SSatish Balay   ilink->nfields = n;
8955a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
8967adad957SLisandro Dalcin   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
8971cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
8985a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
8999005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
90069a612a9SBarry Smith 
901a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
9025a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
9030971522cSBarry Smith 
9040971522cSBarry Smith   if (!next) {
9055a9f2f41SSatish Balay     jac->head       = ilink;
90651f519a2SBarry Smith     ilink->previous = PETSC_NULL;
9070971522cSBarry Smith   } else {
9080971522cSBarry Smith     while (next->next) {
9090971522cSBarry Smith       next = next->next;
9100971522cSBarry Smith     }
9115a9f2f41SSatish Balay     next->next      = ilink;
91251f519a2SBarry Smith     ilink->previous = next;
9130971522cSBarry Smith   }
9140971522cSBarry Smith   jac->nsplits++;
9150971522cSBarry Smith   PetscFunctionReturn(0);
9160971522cSBarry Smith }
9170971522cSBarry Smith EXTERN_C_END
9180971522cSBarry Smith 
919e69d4d44SBarry Smith EXTERN_C_BEGIN
920e69d4d44SBarry Smith #undef __FUNCT__
921e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
9227087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
923e69d4d44SBarry Smith {
924e69d4d44SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
925e69d4d44SBarry Smith   PetscErrorCode ierr;
926e69d4d44SBarry Smith 
927e69d4d44SBarry Smith   PetscFunctionBegin;
928519d70e2SJed Brown   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
929e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
930e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
93113e0d083SBarry Smith   if (n) *n = jac->nsplits;
932e69d4d44SBarry Smith   PetscFunctionReturn(0);
933e69d4d44SBarry Smith }
934e69d4d44SBarry Smith EXTERN_C_END
9350971522cSBarry Smith 
9360971522cSBarry Smith EXTERN_C_BEGIN
9370971522cSBarry Smith #undef __FUNCT__
93869a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
9397087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
9400971522cSBarry Smith {
9410971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
9420971522cSBarry Smith   PetscErrorCode    ierr;
9430971522cSBarry Smith   PetscInt          cnt = 0;
9445a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
9450971522cSBarry Smith 
9460971522cSBarry Smith   PetscFunctionBegin;
9475d480477SMatthew G Knepley   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
9485a9f2f41SSatish Balay   while (ilink) {
9495a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
9505a9f2f41SSatish Balay     ilink = ilink->next;
9510971522cSBarry Smith   }
9525d480477SMatthew 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);
95313e0d083SBarry Smith   if (n) *n = jac->nsplits;
9540971522cSBarry Smith   PetscFunctionReturn(0);
9550971522cSBarry Smith }
9560971522cSBarry Smith EXTERN_C_END
9570971522cSBarry Smith 
958704ba839SBarry Smith EXTERN_C_BEGIN
959704ba839SBarry Smith #undef __FUNCT__
960704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
9617087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
962704ba839SBarry Smith {
963704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
964704ba839SBarry Smith   PetscErrorCode    ierr;
965704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
966704ba839SBarry Smith   char              prefix[128];
967704ba839SBarry Smith 
968704ba839SBarry Smith   PetscFunctionBegin;
9696c924f48SJed Brown   if (jac->splitdefined) {
9706c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
9716c924f48SJed Brown     PetscFunctionReturn(0);
9726c924f48SJed Brown   }
97316913363SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
974a04f6461SBarry Smith   if (splitname) {
975db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
976a04f6461SBarry Smith   } else {
977a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
978a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
979a04f6461SBarry Smith   }
9801ab39975SBarry Smith   ilink->is      = is;
9811ab39975SBarry Smith   ierr           = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
982704ba839SBarry Smith   ilink->next    = PETSC_NULL;
983704ba839SBarry Smith   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
9841cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
985704ba839SBarry Smith   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
9869005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
987704ba839SBarry Smith 
988a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
989704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
990704ba839SBarry Smith 
991704ba839SBarry Smith   if (!next) {
992704ba839SBarry Smith     jac->head       = ilink;
993704ba839SBarry Smith     ilink->previous = PETSC_NULL;
994704ba839SBarry Smith   } else {
995704ba839SBarry Smith     while (next->next) {
996704ba839SBarry Smith       next = next->next;
997704ba839SBarry Smith     }
998704ba839SBarry Smith     next->next      = ilink;
999704ba839SBarry Smith     ilink->previous = next;
1000704ba839SBarry Smith   }
1001704ba839SBarry Smith   jac->nsplits++;
1002704ba839SBarry Smith 
1003704ba839SBarry Smith   PetscFunctionReturn(0);
1004704ba839SBarry Smith }
1005704ba839SBarry Smith EXTERN_C_END
1006704ba839SBarry Smith 
10070971522cSBarry Smith #undef __FUNCT__
10080971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
10090971522cSBarry Smith /*@
10100971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
10110971522cSBarry Smith 
1012ad4df100SBarry Smith     Logically Collective on PC
10130971522cSBarry Smith 
10140971522cSBarry Smith     Input Parameters:
10150971522cSBarry Smith +   pc  - the preconditioner context
1016a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
10170971522cSBarry Smith .   n - the number of fields in this split
1018db4c96c1SJed Brown -   fields - the fields in this split
10190971522cSBarry Smith 
10200971522cSBarry Smith     Level: intermediate
10210971522cSBarry Smith 
1022d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
1023d32f9abdSBarry Smith 
1024d32f9abdSBarry Smith      The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block
1025d32f9abdSBarry 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
1026d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1027d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
1028d32f9abdSBarry Smith 
1029db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
1030db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1031db4c96c1SJed Brown 
1032d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
10330971522cSBarry Smith 
10340971522cSBarry Smith @*/
10357087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields)
10360971522cSBarry Smith {
10374ac538c5SBarry Smith   PetscErrorCode ierr;
10380971522cSBarry Smith 
10390971522cSBarry Smith   PetscFunctionBegin;
10400700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1041db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1042db4c96c1SJed Brown   if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1043db4c96c1SJed Brown   PetscValidIntPointer(fields,3);
10444ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *),(pc,splitname,n,fields));CHKERRQ(ierr);
10450971522cSBarry Smith   PetscFunctionReturn(0);
10460971522cSBarry Smith }
10470971522cSBarry Smith 
10480971522cSBarry Smith #undef __FUNCT__
1049704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
1050704ba839SBarry Smith /*@
1051704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
1052704ba839SBarry Smith 
1053ad4df100SBarry Smith     Logically Collective on PC
1054704ba839SBarry Smith 
1055704ba839SBarry Smith     Input Parameters:
1056704ba839SBarry Smith +   pc  - the preconditioner context
1057a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
1058db4c96c1SJed Brown -   is - the index set that defines the vector elements in this field
1059704ba839SBarry Smith 
1060d32f9abdSBarry Smith 
1061a6ffb8dbSJed Brown     Notes:
1062a6ffb8dbSJed Brown     Use PCFieldSplitSetFields(), for fields defined by strided types.
1063a6ffb8dbSJed Brown 
1064db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
1065db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1066d32f9abdSBarry Smith 
1067704ba839SBarry Smith     Level: intermediate
1068704ba839SBarry Smith 
1069704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
1070704ba839SBarry Smith 
1071704ba839SBarry Smith @*/
10727087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1073704ba839SBarry Smith {
10744ac538c5SBarry Smith   PetscErrorCode ierr;
1075704ba839SBarry Smith 
1076704ba839SBarry Smith   PetscFunctionBegin;
10770700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1078db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1079db4c96c1SJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,3);
10804ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr);
1081704ba839SBarry Smith   PetscFunctionReturn(0);
1082704ba839SBarry Smith }
1083704ba839SBarry Smith 
1084704ba839SBarry Smith #undef __FUNCT__
108557a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS"
108657a9adfeSMatthew G Knepley /*@
108757a9adfeSMatthew G Knepley     PCFieldSplitGetIS - Retrieves the elements for a field as an IS
108857a9adfeSMatthew G Knepley 
108957a9adfeSMatthew G Knepley     Logically Collective on PC
109057a9adfeSMatthew G Knepley 
109157a9adfeSMatthew G Knepley     Input Parameters:
109257a9adfeSMatthew G Knepley +   pc  - the preconditioner context
109357a9adfeSMatthew G Knepley -   splitname - name of this split
109457a9adfeSMatthew G Knepley 
109557a9adfeSMatthew G Knepley     Output Parameter:
109657a9adfeSMatthew G Knepley -   is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found
109757a9adfeSMatthew G Knepley 
109857a9adfeSMatthew G Knepley     Level: intermediate
109957a9adfeSMatthew G Knepley 
110057a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS()
110157a9adfeSMatthew G Knepley 
110257a9adfeSMatthew G Knepley @*/
110357a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
110457a9adfeSMatthew G Knepley {
110557a9adfeSMatthew G Knepley   PetscErrorCode ierr;
110657a9adfeSMatthew G Knepley 
110757a9adfeSMatthew G Knepley   PetscFunctionBegin;
110857a9adfeSMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
110957a9adfeSMatthew G Knepley   PetscValidCharPointer(splitname,2);
111057a9adfeSMatthew G Knepley   PetscValidPointer(is,3);
111157a9adfeSMatthew G Knepley   {
111257a9adfeSMatthew G Knepley     PC_FieldSplit    *jac   = (PC_FieldSplit *) pc->data;
111357a9adfeSMatthew G Knepley     PC_FieldSplitLink ilink = jac->head;
111457a9adfeSMatthew G Knepley     PetscBool         found;
111557a9adfeSMatthew G Knepley 
111657a9adfeSMatthew G Knepley     *is = PETSC_NULL;
111757a9adfeSMatthew G Knepley     while(ilink) {
111857a9adfeSMatthew G Knepley       ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr);
111957a9adfeSMatthew G Knepley       if (found) {
112057a9adfeSMatthew G Knepley         *is = ilink->is;
112157a9adfeSMatthew G Knepley         break;
112257a9adfeSMatthew G Knepley       }
112357a9adfeSMatthew G Knepley       ilink = ilink->next;
112457a9adfeSMatthew G Knepley     }
112557a9adfeSMatthew G Knepley   }
112657a9adfeSMatthew G Knepley   PetscFunctionReturn(0);
112757a9adfeSMatthew G Knepley }
112857a9adfeSMatthew G Knepley 
112957a9adfeSMatthew G Knepley #undef __FUNCT__
113051f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
113151f519a2SBarry Smith /*@
113251f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
113351f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
113451f519a2SBarry Smith 
1135ad4df100SBarry Smith     Logically Collective on PC
113651f519a2SBarry Smith 
113751f519a2SBarry Smith     Input Parameters:
113851f519a2SBarry Smith +   pc  - the preconditioner context
113951f519a2SBarry Smith -   bs - the block size
114051f519a2SBarry Smith 
114151f519a2SBarry Smith     Level: intermediate
114251f519a2SBarry Smith 
114351f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
114451f519a2SBarry Smith 
114551f519a2SBarry Smith @*/
11467087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
114751f519a2SBarry Smith {
11484ac538c5SBarry Smith   PetscErrorCode ierr;
114951f519a2SBarry Smith 
115051f519a2SBarry Smith   PetscFunctionBegin;
11510700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1152c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,bs,2);
11534ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr);
115451f519a2SBarry Smith   PetscFunctionReturn(0);
115551f519a2SBarry Smith }
115651f519a2SBarry Smith 
115751f519a2SBarry Smith #undef __FUNCT__
115869a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
11590971522cSBarry Smith /*@C
116069a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
11610971522cSBarry Smith 
116269a612a9SBarry Smith    Collective on KSP
11630971522cSBarry Smith 
11640971522cSBarry Smith    Input Parameter:
11650971522cSBarry Smith .  pc - the preconditioner context
11660971522cSBarry Smith 
11670971522cSBarry Smith    Output Parameters:
116813e0d083SBarry Smith +  n - the number of splits
116969a612a9SBarry Smith -  pc - the array of KSP contexts
11700971522cSBarry Smith 
11710971522cSBarry Smith    Note:
1172d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
1173d32f9abdSBarry Smith    (not the KSP just the array that contains them).
11740971522cSBarry Smith 
117569a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
11760971522cSBarry Smith 
11770971522cSBarry Smith    Level: advanced
11780971522cSBarry Smith 
11790971522cSBarry Smith .seealso: PCFIELDSPLIT
11800971522cSBarry Smith @*/
11817087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
11820971522cSBarry Smith {
11834ac538c5SBarry Smith   PetscErrorCode ierr;
11840971522cSBarry Smith 
11850971522cSBarry Smith   PetscFunctionBegin;
11860700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
118713e0d083SBarry Smith   if (n) PetscValidIntPointer(n,2);
11884ac538c5SBarry Smith   ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
11890971522cSBarry Smith   PetscFunctionReturn(0);
11900971522cSBarry Smith }
11910971522cSBarry Smith 
1192e69d4d44SBarry Smith #undef __FUNCT__
1193e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition"
1194e69d4d44SBarry Smith /*@
1195e69d4d44SBarry Smith     PCFieldSplitSchurPrecondition -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
1196a04f6461SBarry Smith       A11 matrix. Otherwise no preconditioner is used.
1197e69d4d44SBarry Smith 
1198e69d4d44SBarry Smith     Collective on PC
1199e69d4d44SBarry Smith 
1200e69d4d44SBarry Smith     Input Parameters:
1201e69d4d44SBarry Smith +   pc  - the preconditioner context
1202084e4875SJed Brown .   ptype - which matrix to use for preconditioning the Schur complement
1203084e4875SJed Brown -   userpre - matrix to use for preconditioning, or PETSC_NULL
1204084e4875SJed Brown 
1205084e4875SJed Brown     Notes:
1206a04f6461SBarry Smith     The default is to use the block on the diagonal of the preconditioning matrix.  This is A11, in the (1,1) position.
1207084e4875SJed Brown     There are currently no preconditioners that work directly with the Schur complement so setting
1208084e4875SJed Brown     PC_FIELDSPLIT_SCHUR_PRE_SELF is observationally equivalent to -fieldsplit_1_pc_type none.
1209e69d4d44SBarry Smith 
1210e69d4d44SBarry Smith     Options Database:
1211084e4875SJed Brown .     -pc_fieldsplit_schur_precondition <self,user,diag> default is diag
1212e69d4d44SBarry Smith 
1213e69d4d44SBarry Smith     Level: intermediate
1214e69d4d44SBarry Smith 
1215084e4875SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType
1216e69d4d44SBarry Smith 
1217e69d4d44SBarry Smith @*/
12187087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1219e69d4d44SBarry Smith {
12204ac538c5SBarry Smith   PetscErrorCode ierr;
1221e69d4d44SBarry Smith 
1222e69d4d44SBarry Smith   PetscFunctionBegin;
12230700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12244ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr);
1225e69d4d44SBarry Smith   PetscFunctionReturn(0);
1226e69d4d44SBarry Smith }
1227e69d4d44SBarry Smith 
1228e69d4d44SBarry Smith EXTERN_C_BEGIN
1229e69d4d44SBarry Smith #undef __FUNCT__
1230e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
12317087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1232e69d4d44SBarry Smith {
1233e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1234084e4875SJed Brown   PetscErrorCode  ierr;
1235e69d4d44SBarry Smith 
1236e69d4d44SBarry Smith   PetscFunctionBegin;
1237084e4875SJed Brown   jac->schurpre = ptype;
1238084e4875SJed Brown   if (pre) {
12396bf464f9SBarry Smith     ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
1240084e4875SJed Brown     jac->schur_user = pre;
1241084e4875SJed Brown     ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1242084e4875SJed Brown   }
1243e69d4d44SBarry Smith   PetscFunctionReturn(0);
1244e69d4d44SBarry Smith }
1245e69d4d44SBarry Smith EXTERN_C_END
1246e69d4d44SBarry Smith 
124730ad9308SMatthew Knepley #undef __FUNCT__
124830ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
124930ad9308SMatthew Knepley /*@C
12508c03b21aSDmitry Karpeev    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
125130ad9308SMatthew Knepley 
125230ad9308SMatthew Knepley    Collective on KSP
125330ad9308SMatthew Knepley 
125430ad9308SMatthew Knepley    Input Parameter:
125530ad9308SMatthew Knepley .  pc - the preconditioner context
125630ad9308SMatthew Knepley 
125730ad9308SMatthew Knepley    Output Parameters:
1258a04f6461SBarry Smith +  A00 - the (0,0) block
1259a04f6461SBarry Smith .  A01 - the (0,1) block
1260a04f6461SBarry Smith .  A10 - the (1,0) block
1261a04f6461SBarry Smith -  A11 - the (1,1) block
126230ad9308SMatthew Knepley 
126330ad9308SMatthew Knepley    Level: advanced
126430ad9308SMatthew Knepley 
126530ad9308SMatthew Knepley .seealso: PCFIELDSPLIT
126630ad9308SMatthew Knepley @*/
1267a04f6461SBarry Smith PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
126830ad9308SMatthew Knepley {
126930ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit *) pc->data;
127030ad9308SMatthew Knepley 
127130ad9308SMatthew Knepley   PetscFunctionBegin;
12720700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1273c1235816SBarry Smith   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
1274a04f6461SBarry Smith   if (A00) *A00 = jac->pmat[0];
1275a04f6461SBarry Smith   if (A01) *A01 = jac->B;
1276a04f6461SBarry Smith   if (A10) *A10 = jac->C;
1277a04f6461SBarry Smith   if (A11) *A11 = jac->pmat[1];
127830ad9308SMatthew Knepley   PetscFunctionReturn(0);
127930ad9308SMatthew Knepley }
128030ad9308SMatthew Knepley 
128179416396SBarry Smith EXTERN_C_BEGIN
128279416396SBarry Smith #undef __FUNCT__
128379416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
12847087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
128579416396SBarry Smith {
128679416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1287e69d4d44SBarry Smith   PetscErrorCode ierr;
128879416396SBarry Smith 
128979416396SBarry Smith   PetscFunctionBegin;
129079416396SBarry Smith   jac->type = type;
12913b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
12923b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
12933b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
1294e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1295e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1296e69d4d44SBarry Smith 
12973b224e63SBarry Smith   } else {
12983b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
12993b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
1300e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
13019e7d6b0aSBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr);
13023b224e63SBarry Smith   }
130379416396SBarry Smith   PetscFunctionReturn(0);
130479416396SBarry Smith }
130579416396SBarry Smith EXTERN_C_END
130679416396SBarry Smith 
130751f519a2SBarry Smith EXTERN_C_BEGIN
130851f519a2SBarry Smith #undef __FUNCT__
130951f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
13107087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
131151f519a2SBarry Smith {
131251f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
131351f519a2SBarry Smith 
131451f519a2SBarry Smith   PetscFunctionBegin;
131565e19b50SBarry Smith   if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
131665e19b50SBarry 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);
131751f519a2SBarry Smith   jac->bs = bs;
131851f519a2SBarry Smith   PetscFunctionReturn(0);
131951f519a2SBarry Smith }
132051f519a2SBarry Smith EXTERN_C_END
132151f519a2SBarry Smith 
132279416396SBarry Smith #undef __FUNCT__
132379416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
1324bc08b0f1SBarry Smith /*@
132579416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
132679416396SBarry Smith 
132779416396SBarry Smith    Collective on PC
132879416396SBarry Smith 
132979416396SBarry Smith    Input Parameter:
133079416396SBarry Smith .  pc - the preconditioner context
133181540f2fSBarry Smith .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
133279416396SBarry Smith 
133379416396SBarry Smith    Options Database Key:
1334a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
133579416396SBarry Smith 
133679416396SBarry Smith    Level: Developer
133779416396SBarry Smith 
133879416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
133979416396SBarry Smith 
134079416396SBarry Smith .seealso: PCCompositeSetType()
134179416396SBarry Smith 
134279416396SBarry Smith @*/
13437087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
134479416396SBarry Smith {
13454ac538c5SBarry Smith   PetscErrorCode ierr;
134679416396SBarry Smith 
134779416396SBarry Smith   PetscFunctionBegin;
13480700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13494ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr);
135079416396SBarry Smith   PetscFunctionReturn(0);
135179416396SBarry Smith }
135279416396SBarry Smith 
13530971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
13540971522cSBarry Smith /*MC
1355a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
1356a04f6461SBarry Smith                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.
13570971522cSBarry Smith 
1358edf189efSBarry Smith      To set options on the solvers for each block append -fieldsplit_ to all the PC
1359edf189efSBarry Smith         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
13600971522cSBarry Smith 
1361a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
136269a612a9SBarry Smith          and set the options directly on the resulting KSP object
13630971522cSBarry Smith 
13640971522cSBarry Smith    Level: intermediate
13650971522cSBarry Smith 
136679416396SBarry Smith    Options Database Keys:
136781540f2fSBarry Smith +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
136881540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
136981540f2fSBarry Smith                               been supplied explicitly by -pc_fieldsplit_%d_fields
137081540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
137181540f2fSBarry Smith .   -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative>
1372e69d4d44SBarry Smith .   -pc_fieldsplit_schur_precondition <true,false> default is true
1373435f959eSBarry Smith .   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver
137479416396SBarry Smith 
1375edf189efSBarry Smith -    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
13763b224e63SBarry Smith      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
13773b224e63SBarry Smith 
1378d32f9abdSBarry Smith    Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1379d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
1380d32f9abdSBarry Smith 
1381d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
1382d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1383d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
1384d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
1385d32f9abdSBarry Smith 
1386a04f6461SBarry Smith      For the Schur complement preconditioner if J = ( A00 A01 )
1387a04f6461SBarry Smith                                                     ( A10 A11 )
1388e69d4d44SBarry Smith      the preconditioner is
1389a04f6461SBarry Smith               (I   -B ksp(A00)) ( inv(A00)   0  (I         0  )
1390a04f6461SBarry Smith               (0    I       ) (   0    ksp(S) ) (-A10 ksp(A00) I  )
1391a04f6461SBarry Smith      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of
1392a04f6461SBarry Smith      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).
1393a04f6461SBarry Smith      For PCFieldSplitGetKSP() when field number is
1394edf189efSBarry Smith      0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
1395a04f6461SBarry Smith      A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
13967e8cb189SBarry Smith      option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc
1397e69d4d44SBarry Smith 
1398edf189efSBarry Smith      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
1399edf189efSBarry Smith      is used automatically for a second block.
1400edf189efSBarry Smith 
14010716a85fSBarry Smith      The fieldsplit preconditioner cannot be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1. Generally it should be used with the AIJ format.
14020716a85fSBarry Smith 
1403a541d17aSBarry Smith    Concepts: physics based preconditioners, block preconditioners
14040971522cSBarry Smith 
14057e8cb189SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
1406e69d4d44SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
14070971522cSBarry Smith M*/
14080971522cSBarry Smith 
14090971522cSBarry Smith EXTERN_C_BEGIN
14100971522cSBarry Smith #undef __FUNCT__
14110971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
14127087cfbeSBarry Smith PetscErrorCode  PCCreate_FieldSplit(PC pc)
14130971522cSBarry Smith {
14140971522cSBarry Smith   PetscErrorCode ierr;
14150971522cSBarry Smith   PC_FieldSplit  *jac;
14160971522cSBarry Smith 
14170971522cSBarry Smith   PetscFunctionBegin;
141838f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
14190971522cSBarry Smith   jac->bs        = -1;
14200971522cSBarry Smith   jac->nsplits   = 0;
14213e197d65SBarry Smith   jac->type      = PC_COMPOSITE_MULTIPLICATIVE;
1422e6cab6aaSJed Brown   jac->schurpre  = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
1423c5d2311dSJed Brown   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL;
142451f519a2SBarry Smith 
14250971522cSBarry Smith   pc->data     = (void*)jac;
14260971522cSBarry Smith 
14270971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
1428421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
14290971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
1430574deadeSBarry Smith   pc->ops->reset             = PCReset_FieldSplit;
14310971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
14320971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
14330971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
14340971522cSBarry Smith   pc->ops->applyrichardson   = 0;
14350971522cSBarry Smith 
143669a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
143769a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
14380971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
14390971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1440704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
1441704ba839SBarry Smith                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
144279416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
144379416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
144451f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
144551f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
14460971522cSBarry Smith   PetscFunctionReturn(0);
14470971522cSBarry Smith }
14480971522cSBarry Smith EXTERN_C_END
14490971522cSBarry Smith 
14500971522cSBarry Smith 
1451a541d17aSBarry Smith 
1452