xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision ea6813d4b88b8c416dfff52ac933cea8baff5de8)
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;
4863ec74ffSBarry Smith   PetscBool                          reset;         /* indicates PCReset() has been last called on this object, hack */
49c1570756SJed Brown   PetscBool                          suboptionsset; /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */
500971522cSBarry Smith } PC_FieldSplit;
510971522cSBarry Smith 
5216913363SBarry Smith /*
5316913363SBarry Smith     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
5416913363SBarry Smith    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
5516913363SBarry Smith    PC you could change this.
5616913363SBarry Smith */
57084e4875SJed Brown 
58e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it.  This way the
59084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */
60084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac)
61084e4875SJed Brown {
62084e4875SJed Brown   switch (jac->schurpre) {
63084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur;
64084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1];
65084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */
66084e4875SJed Brown     default:
67084e4875SJed Brown       return jac->schur_user ? jac->schur_user : jac->pmat[1];
68084e4875SJed Brown   }
69084e4875SJed Brown }
70084e4875SJed Brown 
71084e4875SJed Brown 
720971522cSBarry Smith #undef __FUNCT__
730971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
740971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
750971522cSBarry Smith {
760971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
770971522cSBarry Smith   PetscErrorCode    ierr;
78ace3abfcSBarry Smith   PetscBool         iascii;
790971522cSBarry Smith   PetscInt          i,j;
805a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
810971522cSBarry Smith 
820971522cSBarry Smith   PetscFunctionBegin;
832692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
840971522cSBarry Smith   if (iascii) {
8551f519a2SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
8669a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
870971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
880971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
891ab39975SBarry Smith       if (ilink->fields) {
900971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
9179416396SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
925a9f2f41SSatish Balay 	for (j=0; j<ilink->nfields; j++) {
9379416396SBarry Smith 	  if (j > 0) {
9479416396SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
9579416396SBarry Smith 	  }
965a9f2f41SSatish Balay 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
970971522cSBarry Smith 	}
980971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
9979416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1001ab39975SBarry Smith       } else {
1011ab39975SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1021ab39975SBarry Smith       }
1035a9f2f41SSatish Balay       ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
1045a9f2f41SSatish Balay       ilink = ilink->next;
1050971522cSBarry Smith     }
1060971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1070971522cSBarry Smith   } else {
10865e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1090971522cSBarry Smith   }
1100971522cSBarry Smith   PetscFunctionReturn(0);
1110971522cSBarry Smith }
1120971522cSBarry Smith 
1130971522cSBarry Smith #undef __FUNCT__
1143b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur"
1153b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
1163b224e63SBarry Smith {
1173b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1183b224e63SBarry Smith   PetscErrorCode    ierr;
119ace3abfcSBarry Smith   PetscBool         iascii;
1203b224e63SBarry Smith   PetscInt          i,j;
1213b224e63SBarry Smith   PC_FieldSplitLink ilink = jac->head;
1223b224e63SBarry Smith   KSP               ksp;
1233b224e63SBarry Smith 
1243b224e63SBarry Smith   PetscFunctionBegin;
1252692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1263b224e63SBarry Smith   if (iascii) {
127c5d2311dSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactorizationTypes[jac->schurfactorization]);CHKERRQ(ierr);
1283b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Split info:\n");CHKERRQ(ierr);
1293b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1303b224e63SBarry Smith     for (i=0; i<jac->nsplits; i++) {
1313b224e63SBarry Smith       if (ilink->fields) {
1323b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
1333b224e63SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
1343b224e63SBarry Smith 	for (j=0; j<ilink->nfields; j++) {
1353b224e63SBarry Smith 	  if (j > 0) {
1363b224e63SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
1373b224e63SBarry Smith 	  }
1383b224e63SBarry Smith 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
1393b224e63SBarry Smith 	}
1403b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1413b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1423b224e63SBarry Smith       } else {
1433b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1443b224e63SBarry Smith       }
1453b224e63SBarry Smith       ilink = ilink->next;
1463b224e63SBarry Smith     }
147435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block \n");CHKERRQ(ierr);
1483b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
14912cae6f2SJed Brown     if (jac->schur) {
15012cae6f2SJed Brown       ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
1513b224e63SBarry Smith       ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
15212cae6f2SJed Brown     } else {
15312cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
15412cae6f2SJed Brown     }
1553b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
156435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr);
1573b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15812cae6f2SJed Brown     if (jac->kspschur) {
1593b224e63SBarry Smith       ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
16012cae6f2SJed Brown     } else {
16112cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
16212cae6f2SJed Brown     }
1633b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1643b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1653b224e63SBarry Smith   } else {
16665e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1673b224e63SBarry Smith   }
1683b224e63SBarry Smith   PetscFunctionReturn(0);
1693b224e63SBarry Smith }
1703b224e63SBarry Smith 
1713b224e63SBarry Smith #undef __FUNCT__
1726c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private"
1736c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */
1746c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc)
1756c924f48SJed Brown {
1766c924f48SJed Brown   PetscErrorCode ierr;
1776c924f48SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1786c924f48SJed Brown   PetscInt       i,nfields,*ifields;
179ace3abfcSBarry Smith   PetscBool      flg;
1806c924f48SJed Brown   char           optionname[128],splitname[8];
1816c924f48SJed Brown 
1826c924f48SJed Brown   PetscFunctionBegin;
1836c924f48SJed Brown   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr);
1846c924f48SJed Brown   for (i=0,flg=PETSC_TRUE; ; i++) {
1856c924f48SJed Brown     ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
1866c924f48SJed Brown     ierr = PetscSNPrintf(optionname,sizeof optionname,"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr);
1876c924f48SJed Brown     nfields = jac->bs;
1886c924f48SJed Brown     ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr);
1896c924f48SJed Brown     if (!flg) break;
190*ea6813d4SMatthew G Knepley     if (!nfields) SETERRQ(((PetscObject) pc)->comm,PETSC_ERR_USER,"Cannot list zero fields");
1916c924f48SJed Brown     ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields);CHKERRQ(ierr);
1926c924f48SJed Brown   }
1936c924f48SJed Brown   if (i > 0) {
1946c924f48SJed Brown     /* Makes command-line setting of splits take precedence over setting them in code.
1956c924f48SJed Brown        Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would
1966c924f48SJed Brown        create new splits, which would probably not be what the user wanted. */
1976c924f48SJed Brown     jac->splitdefined = PETSC_TRUE;
1986c924f48SJed Brown   }
1996c924f48SJed Brown   ierr = PetscFree(ifields);CHKERRQ(ierr);
2006c924f48SJed Brown   PetscFunctionReturn(0);
2016c924f48SJed Brown }
2026c924f48SJed Brown 
2036c924f48SJed Brown #undef __FUNCT__
20469a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
20569a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
2060971522cSBarry Smith {
2070971522cSBarry Smith   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
2080971522cSBarry Smith   PetscErrorCode    ierr;
2095a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
2106ce1633cSBarry Smith   PetscBool         flg = PETSC_FALSE,stokes = PETSC_FALSE;
2116c924f48SJed Brown   PetscInt          i;
2120971522cSBarry Smith 
2130971522cSBarry Smith   PetscFunctionBegin;
214d32f9abdSBarry Smith   if (!ilink) {
215d0af7cd3SBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
216d0af7cd3SBarry Smith     if (pc->dm && !stokes) {
217*ea6813d4SMatthew G Knepley       PetscBool dmcomposite, dmcomplex;
2188b8307b2SJed Brown       ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr);
219*ea6813d4SMatthew G Knepley       ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPLEX,&dmcomplex);CHKERRQ(ierr);
2208b8307b2SJed Brown       if (dmcomposite) {
2218b8307b2SJed Brown         PetscInt nDM;
2228b8307b2SJed Brown         IS       *fields;
2232fa5ba8aSJed Brown         DM       *dms;
2248b8307b2SJed Brown         ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
2258b8307b2SJed Brown         ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr);
2268b8307b2SJed Brown         ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr);
2272fa5ba8aSJed Brown         ierr = PetscMalloc(nDM*sizeof(DM),&dms);CHKERRQ(ierr);
2282fa5ba8aSJed Brown         ierr = DMCompositeGetEntriesArray(pc->dm,dms);CHKERRQ(ierr);
2298b8307b2SJed Brown         for (i=0; i<nDM; i++) {
2302fa5ba8aSJed Brown           char buf[256];
2312fa5ba8aSJed Brown           const char *splitname;
2322fa5ba8aSJed Brown           /* Split naming precedence: object name, prefix, number */
2332fa5ba8aSJed Brown           splitname = ((PetscObject)pc->dm)->name;
2342fa5ba8aSJed Brown           if (!splitname) {
2352fa5ba8aSJed Brown             ierr = PetscObjectGetOptionsPrefix((PetscObject)dms[i],&splitname);CHKERRQ(ierr);
2362fa5ba8aSJed Brown             if (splitname) {
2372fa5ba8aSJed Brown               size_t len;
2382fa5ba8aSJed Brown               ierr = PetscStrncpy(buf,splitname,sizeof buf);CHKERRQ(ierr);
2392fa5ba8aSJed Brown               buf[sizeof buf - 1] = 0;
2402fa5ba8aSJed Brown               ierr = PetscStrlen(buf,&len);CHKERRQ(ierr);
2412fa5ba8aSJed Brown               if (buf[len-1] == '_') buf[len-1] = 0; /* Remove trailing underscore if it was used */
2422fa5ba8aSJed Brown               splitname = buf;
2432fa5ba8aSJed Brown             }
2442fa5ba8aSJed Brown           }
2452fa5ba8aSJed Brown           if (!splitname) {
2462fa5ba8aSJed Brown             ierr = PetscSNPrintf(buf,sizeof buf,"%D",i);CHKERRQ(ierr);
2472fa5ba8aSJed Brown             splitname = buf;
2482fa5ba8aSJed Brown           }
2498b8307b2SJed Brown           ierr = PCFieldSplitSetIS(pc,splitname,fields[i]);CHKERRQ(ierr);
250fcfd50ebSBarry Smith           ierr = ISDestroy(&fields[i]);CHKERRQ(ierr);
2518b8307b2SJed Brown         }
2528b8307b2SJed Brown         ierr = PetscFree(fields);CHKERRQ(ierr);
2532fa5ba8aSJed Brown         for (ilink=jac->head,i=0; ilink; ilink=ilink->next,i++) {
2542fa5ba8aSJed Brown           ierr = KSPSetDM(ilink->ksp,dms[i]);CHKERRQ(ierr);
2552fa5ba8aSJed Brown           ierr = KSPSetDMActive(ilink->ksp,PETSC_FALSE);CHKERRQ(ierr);
2562fa5ba8aSJed Brown         }
2572fa5ba8aSJed Brown         ierr = PetscFree(dms);CHKERRQ(ierr);
258*ea6813d4SMatthew G Knepley       } else if (dmcomplex) {
259*ea6813d4SMatthew G Knepley         /* Define fields */
260*ea6813d4SMatthew G Knepley         SETERRQ(((PetscObject) pc)->comm, PETSC_ERR_SUP, "IMPLEMENT THIS");
2618b8307b2SJed Brown       }
26266ffff09SJed Brown     } else {
263521d7252SBarry Smith       if (jac->bs <= 0) {
264704ba839SBarry Smith         if (pc->pmat) {
265521d7252SBarry Smith           ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
266704ba839SBarry Smith         } else {
267704ba839SBarry Smith           jac->bs = 1;
268704ba839SBarry Smith         }
269521d7252SBarry Smith       }
270d32f9abdSBarry Smith 
271acfcf0e5SJed Brown       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
2726ce1633cSBarry Smith       if (stokes) {
2736ce1633cSBarry Smith         IS       zerodiags,rest;
2746ce1633cSBarry Smith         PetscInt nmin,nmax;
2756ce1633cSBarry Smith 
2766ce1633cSBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
2776ce1633cSBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
2786ce1633cSBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
2796ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
2806ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr);
281fcfd50ebSBarry Smith         ierr = ISDestroy(&zerodiags);CHKERRQ(ierr);
282fcfd50ebSBarry Smith         ierr = ISDestroy(&rest);CHKERRQ(ierr);
2836ce1633cSBarry Smith       } else {
284d32f9abdSBarry Smith         if (!flg) {
285d32f9abdSBarry Smith           /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
286d32f9abdSBarry Smith            then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
2876c924f48SJed Brown           ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
2886c924f48SJed Brown           if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
289d32f9abdSBarry Smith         }
2906c924f48SJed Brown         if (flg || !jac->splitdefined) {
291d32f9abdSBarry Smith           ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
292db4c96c1SJed Brown           for (i=0; i<jac->bs; i++) {
2936c924f48SJed Brown             char splitname[8];
2946c924f48SJed Brown             ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
295db4c96c1SJed Brown             ierr = PCFieldSplitSetFields(pc,splitname,1,&i);CHKERRQ(ierr);
29679416396SBarry Smith           }
29797bbdb24SBarry Smith           jac->defaultsplit = PETSC_TRUE;
298521d7252SBarry Smith         }
29966ffff09SJed Brown       }
3006ce1633cSBarry Smith     }
301edf189efSBarry Smith   } else if (jac->nsplits == 1) {
302edf189efSBarry Smith     if (ilink->is) {
303edf189efSBarry Smith       IS       is2;
304edf189efSBarry Smith       PetscInt nmin,nmax;
305edf189efSBarry Smith 
306edf189efSBarry Smith       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
307edf189efSBarry Smith       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
308db4c96c1SJed Brown       ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr);
309fcfd50ebSBarry Smith       ierr = ISDestroy(&is2);CHKERRQ(ierr);
310*ea6813d4SMatthew G Knepley     } else SETERRQ(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
31163ec74ffSBarry Smith   } else if (jac->reset) {
31263ec74ffSBarry Smith     /* PCReset() has been called on this PC, ilink exists but all IS and Vec data structures in it must be rebuilt
313d0af7cd3SBarry Smith        This is basically the !ilink portion of code above copied from above and the allocation of the ilinks removed
314d0af7cd3SBarry Smith        since they already exist. This should be totally rewritten */
315d0af7cd3SBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
316d0af7cd3SBarry Smith     if (pc->dm && !stokes) {
317d0af7cd3SBarry Smith       PetscBool dmcomposite;
318d0af7cd3SBarry Smith       ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr);
319d0af7cd3SBarry Smith       if (dmcomposite) {
320d0af7cd3SBarry Smith         PetscInt nDM;
321d0af7cd3SBarry Smith         IS       *fields;
3222fa5ba8aSJed Brown         DM       *dms;
323d0af7cd3SBarry Smith         ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
324d0af7cd3SBarry Smith         ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr);
325d0af7cd3SBarry Smith         ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr);
3262fa5ba8aSJed Brown         ierr = PetscMalloc(nDM*sizeof(DM),&dms);CHKERRQ(ierr);
3272fa5ba8aSJed Brown         ierr = DMCompositeGetEntriesArray(pc->dm,dms);CHKERRQ(ierr);
328d0af7cd3SBarry Smith         for (i=0; i<nDM; i++) {
3292fa5ba8aSJed Brown           ierr = KSPSetDM(ilink->ksp,dms[i]);CHKERRQ(ierr);
3302fa5ba8aSJed Brown           ierr = KSPSetDMActive(ilink->ksp,PETSC_FALSE);CHKERRQ(ierr);
331d0af7cd3SBarry Smith           ilink->is = fields[i];
332d0af7cd3SBarry Smith           ilink     = ilink->next;
333edf189efSBarry Smith         }
334d0af7cd3SBarry Smith         ierr = PetscFree(fields);CHKERRQ(ierr);
3352fa5ba8aSJed Brown         ierr = PetscFree(dms);CHKERRQ(ierr);
336d0af7cd3SBarry Smith       }
337d0af7cd3SBarry Smith     } else {
338d0af7cd3SBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
339d0af7cd3SBarry Smith       if (stokes) {
340d0af7cd3SBarry Smith         IS       zerodiags,rest;
341d0af7cd3SBarry Smith         PetscInt nmin,nmax;
342d0af7cd3SBarry Smith 
343d0af7cd3SBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
344d0af7cd3SBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
345d0af7cd3SBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
34620b26d62SBarry Smith         ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
34720b26d62SBarry Smith         ierr = ISDestroy(&ilink->next->is);CHKERRQ(ierr);
348d0af7cd3SBarry Smith         ilink->is       = rest;
349d0af7cd3SBarry Smith         ilink->next->is = zerodiags;
35020b26d62SBarry Smith       } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used");
351d0af7cd3SBarry Smith     }
352d0af7cd3SBarry Smith   }
353d0af7cd3SBarry Smith 
354*ea6813d4SMatthew G Knepley   if (jac->nsplits < 2) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits);
35569a612a9SBarry Smith   PetscFunctionReturn(0);
35669a612a9SBarry Smith }
35769a612a9SBarry Smith 
35869a612a9SBarry Smith #undef __FUNCT__
35969a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
36069a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
36169a612a9SBarry Smith {
36269a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
36369a612a9SBarry Smith   PetscErrorCode    ierr;
3645a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
365cf502942SBarry Smith   PetscInt          i,nsplit,ccsize;
36669a612a9SBarry Smith   MatStructure      flag = pc->flag;
367ace3abfcSBarry Smith   PetscBool         sorted;
36869a612a9SBarry Smith 
36969a612a9SBarry Smith   PetscFunctionBegin;
37069a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
37197bbdb24SBarry Smith   nsplit = jac->nsplits;
3725a9f2f41SSatish Balay   ilink  = jac->head;
37397bbdb24SBarry Smith 
37497bbdb24SBarry Smith   /* get the matrices for each split */
375704ba839SBarry Smith   if (!jac->issetup) {
3761b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
37797bbdb24SBarry Smith 
378704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
379704ba839SBarry Smith 
380704ba839SBarry Smith     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
38151f519a2SBarry Smith     bs     = jac->bs;
38297bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
383cf502942SBarry Smith     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
3841b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
3851b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
3861b9fc7fcSBarry Smith       if (jac->defaultsplit) {
387704ba839SBarry Smith         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
388704ba839SBarry Smith       } else if (!ilink->is) {
389ccb205f8SBarry Smith         if (ilink->nfields > 1) {
3905a9f2f41SSatish Balay           PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
3915a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
3921b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
3931b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
3941b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
39597bbdb24SBarry Smith             }
39697bbdb24SBarry Smith           }
397d67e408aSBarry Smith           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr);
398ccb205f8SBarry Smith         } else {
399704ba839SBarry Smith           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
400ccb205f8SBarry Smith         }
4013e197d65SBarry Smith       }
402704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
403e32f2f54SBarry Smith       if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
404704ba839SBarry Smith       ilink = ilink->next;
4051b9fc7fcSBarry Smith     }
4061b9fc7fcSBarry Smith   }
4071b9fc7fcSBarry Smith 
408704ba839SBarry Smith   ilink  = jac->head;
40997bbdb24SBarry Smith   if (!jac->pmat) {
410cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
411cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
4124aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
413704ba839SBarry Smith       ilink = ilink->next;
414cf502942SBarry Smith     }
41597bbdb24SBarry Smith   } else {
416cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
4174aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
418704ba839SBarry Smith       ilink = ilink->next;
419cf502942SBarry Smith     }
42097bbdb24SBarry Smith   }
421519d70e2SJed Brown   if (jac->realdiagonal) {
422519d70e2SJed Brown     ilink = jac->head;
423519d70e2SJed Brown     if (!jac->mat) {
424519d70e2SJed Brown       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr);
425519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
426519d70e2SJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
427519d70e2SJed Brown         ilink = ilink->next;
428519d70e2SJed Brown       }
429519d70e2SJed Brown     } else {
430519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
431966e74d7SJed Brown         if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);}
432519d70e2SJed Brown         ilink = ilink->next;
433519d70e2SJed Brown       }
434519d70e2SJed Brown     }
435519d70e2SJed Brown   } else {
436519d70e2SJed Brown     jac->mat = jac->pmat;
437519d70e2SJed Brown   }
43897bbdb24SBarry Smith 
4396c8605c2SJed Brown   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
44068dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
44168dd23aaSBarry Smith     ilink  = jac->head;
44268dd23aaSBarry Smith     if (!jac->Afield) {
44368dd23aaSBarry Smith       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
44468dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
4454aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
44668dd23aaSBarry Smith         ilink = ilink->next;
44768dd23aaSBarry Smith       }
44868dd23aaSBarry Smith     } else {
44968dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
4504aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
45168dd23aaSBarry Smith         ilink = ilink->next;
45268dd23aaSBarry Smith       }
45368dd23aaSBarry Smith     }
45468dd23aaSBarry Smith   }
45568dd23aaSBarry Smith 
4563b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
4573b224e63SBarry Smith     IS       ccis;
4584aa3045dSJed Brown     PetscInt rstart,rend;
459e7e72b3dSBarry Smith     if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
46068dd23aaSBarry Smith 
461e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
462e6cab6aaSJed Brown     ierr  = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
463e6cab6aaSJed Brown 
4643b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
4653b224e63SBarry Smith     if (jac->schur) {
4663b224e63SBarry Smith       ilink = jac->head;
4674aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4684aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
469fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4703b224e63SBarry Smith       ilink = ilink->next;
4714aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4724aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
473fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
474519d70e2SJed Brown       ierr  = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr);
475084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr);
4763b224e63SBarry Smith 
4773b224e63SBarry Smith      } else {
4781cee3971SBarry Smith       KSP ksp;
4796c924f48SJed Brown       char schurprefix[256];
4803b224e63SBarry Smith 
481a04f6461SBarry Smith       /* extract the A01 and A10 matrices */
4823b224e63SBarry Smith       ilink = jac->head;
4834aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4844aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
485fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4863b224e63SBarry Smith       ilink = ilink->next;
4874aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4884aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
489fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4907d50072eSJed Brown       /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */
491519d70e2SJed Brown       ierr  = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr);
492a04f6461SBarry Smith       /* set tabbing and options prefix of KSP inside the MatSchur */
4931cee3971SBarry Smith       ierr  = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
494e69d4d44SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr);
495a04f6461SBarry Smith       ierr  = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",jac->head->splitname);CHKERRQ(ierr);
496a04f6461SBarry Smith       ierr  = KSPSetOptionsPrefix(ksp,schurprefix);CHKERRQ(ierr);
49720b26d62SBarry Smith       /* Need to call this everytime because new matrix is being created */
49820b26d62SBarry Smith       ierr  = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
4997ae8954aSSatish Balay       ierr  = MatSetUp(jac->schur);CHKERRQ(ierr);
5003b224e63SBarry Smith 
5013b224e63SBarry Smith       ierr  = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr);
5029005cf84SBarry Smith       ierr  = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr);
5031cee3971SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
504084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
505084e4875SJed Brown       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
506084e4875SJed Brown         PC pc;
507cd5f4a64SJed Brown         ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr);
508084e4875SJed Brown         ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
509084e4875SJed Brown         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
510e69d4d44SBarry Smith       }
5116c924f48SJed Brown       ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
5126c924f48SJed Brown       ierr  = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr);
5133b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
51420b26d62SBarry Smith       /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */
51520b26d62SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
5163b224e63SBarry Smith 
5173b224e63SBarry Smith       ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr);
5183b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr);
5193b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr);
5203b224e63SBarry Smith       ilink = jac->head;
5213b224e63SBarry Smith       ilink->x = jac->x[0]; ilink->y = jac->y[0];
5223b224e63SBarry Smith       ilink = ilink->next;
5233b224e63SBarry Smith       ilink->x = jac->x[1]; ilink->y = jac->y[1];
5243b224e63SBarry Smith     }
5253b224e63SBarry Smith   } else {
52697bbdb24SBarry Smith     /* set up the individual PCs */
52797bbdb24SBarry Smith     i    = 0;
5285a9f2f41SSatish Balay     ilink = jac->head;
5295a9f2f41SSatish Balay     while (ilink) {
530519d70e2SJed Brown       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr);
5313b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
532c1570756SJed Brown       if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);}
5335a9f2f41SSatish Balay       ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
53497bbdb24SBarry Smith       i++;
5355a9f2f41SSatish Balay       ilink = ilink->next;
5360971522cSBarry Smith     }
53797bbdb24SBarry Smith 
53897bbdb24SBarry Smith     /* create work vectors for each split */
5391b9fc7fcSBarry Smith     if (!jac->x) {
54097bbdb24SBarry Smith       ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
5415a9f2f41SSatish Balay       ilink = jac->head;
54297bbdb24SBarry Smith       for (i=0; i<nsplit; i++) {
543906ed7ccSBarry Smith         Vec *vl,*vr;
544906ed7ccSBarry Smith 
545906ed7ccSBarry Smith         ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
546906ed7ccSBarry Smith         ilink->x  = *vr;
547906ed7ccSBarry Smith         ilink->y  = *vl;
548906ed7ccSBarry Smith         ierr      = PetscFree(vr);CHKERRQ(ierr);
549906ed7ccSBarry Smith         ierr      = PetscFree(vl);CHKERRQ(ierr);
5505a9f2f41SSatish Balay         jac->x[i] = ilink->x;
5515a9f2f41SSatish Balay         jac->y[i] = ilink->y;
5525a9f2f41SSatish Balay         ilink     = ilink->next;
55397bbdb24SBarry Smith       }
5543b224e63SBarry Smith     }
5553b224e63SBarry Smith   }
5563b224e63SBarry Smith 
5573b224e63SBarry Smith 
558d0f46423SBarry Smith   if (!jac->head->sctx) {
5593b224e63SBarry Smith     Vec xtmp;
5603b224e63SBarry Smith 
56179416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
5621b9fc7fcSBarry Smith 
5635a9f2f41SSatish Balay     ilink = jac->head;
5641b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
5651b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
566704ba839SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
5675a9f2f41SSatish Balay       ilink = ilink->next;
5681b9fc7fcSBarry Smith     }
5696bf464f9SBarry Smith     ierr = VecDestroy(&xtmp);CHKERRQ(ierr);
5701b9fc7fcSBarry Smith   }
571c1570756SJed Brown   jac->suboptionsset = PETSC_TRUE;
5720971522cSBarry Smith   PetscFunctionReturn(0);
5730971522cSBarry Smith }
5740971522cSBarry Smith 
5755a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
576ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
577ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
5785a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
579ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
580ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
58179416396SBarry Smith 
5820971522cSBarry Smith #undef __FUNCT__
5833b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur"
5843b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
5853b224e63SBarry Smith {
5863b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
5873b224e63SBarry Smith   PetscErrorCode    ierr;
5883b224e63SBarry Smith   KSP               ksp;
5893b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
5903b224e63SBarry Smith 
5913b224e63SBarry Smith   PetscFunctionBegin;
5923b224e63SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
5933b224e63SBarry Smith 
594c5d2311dSJed Brown   switch (jac->schurfactorization) {
595c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG:
596a04f6461SBarry Smith       /* [A00 0; 0 -S], positive definite, suitable for MINRES */
597c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
598c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
599c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
600c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
601c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
602c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
603c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
604c5d2311dSJed Brown       ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr);
605c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
606c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
607c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
608c5d2311dSJed Brown       break;
609c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER:
610a04f6461SBarry Smith       /* [A00 0; A10 S], suitable for left preconditioning */
611c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
612c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
613c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
614c5d2311dSJed Brown       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
615c5d2311dSJed Brown       ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr);
616c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
617c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
618c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
619c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
620c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
621c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
622c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
623c5d2311dSJed Brown       break;
624c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER:
625a04f6461SBarry Smith       /* [A00 A01; 0 S], suitable for right preconditioning */
626c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
627c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
628c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
629c5d2311dSJed Brown       ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
630c5d2311dSJed Brown       ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr);
631c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
632c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
633c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
634c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
635c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
636c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
637c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
638c5d2311dSJed Brown       break;
639c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL:
640a04f6461SBarry 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 */
6413b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6423b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6433b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
6443b224e63SBarry Smith       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
6453b224e63SBarry Smith       ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
6463b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6473b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6483b224e63SBarry Smith 
6493b224e63SBarry Smith       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
6503b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6513b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6523b224e63SBarry Smith 
6533b224e63SBarry Smith       ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
6543b224e63SBarry Smith       ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
6553b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
6563b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6573b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
658c5d2311dSJed Brown   }
6593b224e63SBarry Smith   PetscFunctionReturn(0);
6603b224e63SBarry Smith }
6613b224e63SBarry Smith 
6623b224e63SBarry Smith #undef __FUNCT__
6630971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
6640971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
6650971522cSBarry Smith {
6660971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
6670971522cSBarry Smith   PetscErrorCode    ierr;
6685a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
669af93645bSJed Brown   PetscInt          cnt;
6700971522cSBarry Smith 
6710971522cSBarry Smith   PetscFunctionBegin;
67251f519a2SBarry Smith   CHKMEMQ;
67351f519a2SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
67451f519a2SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
67551f519a2SBarry Smith 
67679416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
6771b9fc7fcSBarry Smith     if (jac->defaultsplit) {
6780971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
6795a9f2f41SSatish Balay       while (ilink) {
6805a9f2f41SSatish Balay         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
6815a9f2f41SSatish Balay         ilink = ilink->next;
6820971522cSBarry Smith       }
6830971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
6841b9fc7fcSBarry Smith     } else {
685efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
6865a9f2f41SSatish Balay       while (ilink) {
6875a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
6885a9f2f41SSatish Balay         ilink = ilink->next;
6891b9fc7fcSBarry Smith       }
6901b9fc7fcSBarry Smith     }
69116913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
69279416396SBarry Smith     if (!jac->w1) {
69379416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
69479416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
69579416396SBarry Smith     }
696efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
6975a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
6983e197d65SBarry Smith     cnt = 1;
6995a9f2f41SSatish Balay     while (ilink->next) {
7005a9f2f41SSatish Balay       ilink = ilink->next;
7013e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
7023e197d65SBarry Smith       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
7033e197d65SBarry Smith       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
7043e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7053e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7063e197d65SBarry Smith       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
7073e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7083e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7093e197d65SBarry Smith     }
71051f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
71111755939SBarry Smith       cnt -= 2;
71251f519a2SBarry Smith       while (ilink->previous) {
71351f519a2SBarry Smith         ilink = ilink->previous;
71411755939SBarry Smith         /* compute the residual only over the part of the vector needed */
71511755939SBarry Smith         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
71611755939SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
71711755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
71811755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
71911755939SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
72011755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
72111755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
72251f519a2SBarry Smith       }
72311755939SBarry Smith     }
72465e19b50SBarry Smith   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
72551f519a2SBarry Smith   CHKMEMQ;
7260971522cSBarry Smith   PetscFunctionReturn(0);
7270971522cSBarry Smith }
7280971522cSBarry Smith 
729421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
730ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
731ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
732421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
733ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
734ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
735421e10b8SBarry Smith 
736421e10b8SBarry Smith #undef __FUNCT__
7378c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit"
738421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
739421e10b8SBarry Smith {
740421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
741421e10b8SBarry Smith   PetscErrorCode    ierr;
742421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
743421e10b8SBarry Smith 
744421e10b8SBarry Smith   PetscFunctionBegin;
745421e10b8SBarry Smith   CHKMEMQ;
746421e10b8SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
747421e10b8SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
748421e10b8SBarry Smith 
749421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
750421e10b8SBarry Smith     if (jac->defaultsplit) {
751421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
752421e10b8SBarry Smith       while (ilink) {
753421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
754421e10b8SBarry Smith 	ilink = ilink->next;
755421e10b8SBarry Smith       }
756421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
757421e10b8SBarry Smith     } else {
758421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
759421e10b8SBarry Smith       while (ilink) {
760421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
761421e10b8SBarry Smith 	ilink = ilink->next;
762421e10b8SBarry Smith       }
763421e10b8SBarry Smith     }
764421e10b8SBarry Smith   } else {
765421e10b8SBarry Smith     if (!jac->w1) {
766421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
767421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
768421e10b8SBarry Smith     }
769421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
770421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
771421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
772421e10b8SBarry Smith       while (ilink->next) {
773421e10b8SBarry Smith         ilink = ilink->next;
7749989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
775421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
776421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
777421e10b8SBarry Smith       }
778421e10b8SBarry Smith       while (ilink->previous) {
779421e10b8SBarry Smith         ilink = ilink->previous;
7809989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
781421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
782421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
783421e10b8SBarry Smith       }
784421e10b8SBarry Smith     } else {
785421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
786421e10b8SBarry Smith 	ilink = ilink->next;
787421e10b8SBarry Smith       }
788421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
789421e10b8SBarry Smith       while (ilink->previous) {
790421e10b8SBarry Smith 	ilink = ilink->previous;
7919989ab13SBarry Smith 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
792421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
793421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
794421e10b8SBarry Smith       }
795421e10b8SBarry Smith     }
796421e10b8SBarry Smith   }
797421e10b8SBarry Smith   CHKMEMQ;
798421e10b8SBarry Smith   PetscFunctionReturn(0);
799421e10b8SBarry Smith }
800421e10b8SBarry Smith 
8010971522cSBarry Smith #undef __FUNCT__
802574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit"
803574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc)
8040971522cSBarry Smith {
8050971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
8060971522cSBarry Smith   PetscErrorCode    ierr;
8075a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
8080971522cSBarry Smith 
8090971522cSBarry Smith   PetscFunctionBegin;
8105a9f2f41SSatish Balay   while (ilink) {
811574deadeSBarry Smith     ierr = KSPReset(ilink->ksp);CHKERRQ(ierr);
812fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->x);CHKERRQ(ierr);
813fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->y);CHKERRQ(ierr);
814fcfd50ebSBarry Smith     ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr);
815fcfd50ebSBarry Smith     ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
8165a9f2f41SSatish Balay     next = ilink->next;
8175a9f2f41SSatish Balay     ilink = next;
8180971522cSBarry Smith   }
81905b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
820574deadeSBarry Smith   if (jac->mat && jac->mat != jac->pmat) {
821574deadeSBarry Smith     ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);
822574deadeSBarry Smith   } else if (jac->mat) {
823574deadeSBarry Smith     jac->mat = PETSC_NULL;
824574deadeSBarry Smith   }
82597bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
82668dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
8276bf464f9SBarry Smith   ierr = VecDestroy(&jac->w1);CHKERRQ(ierr);
8286bf464f9SBarry Smith   ierr = VecDestroy(&jac->w2);CHKERRQ(ierr);
8296bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur);CHKERRQ(ierr);
8306bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
8316bf464f9SBarry Smith   ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr);
8326bf464f9SBarry Smith   ierr = MatDestroy(&jac->B);CHKERRQ(ierr);
8336bf464f9SBarry Smith   ierr = MatDestroy(&jac->C);CHKERRQ(ierr);
83463ec74ffSBarry Smith   jac->reset = PETSC_TRUE;
835574deadeSBarry Smith   PetscFunctionReturn(0);
836574deadeSBarry Smith }
837574deadeSBarry Smith 
838574deadeSBarry Smith #undef __FUNCT__
839574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
840574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
841574deadeSBarry Smith {
842574deadeSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
843574deadeSBarry Smith   PetscErrorCode    ierr;
844574deadeSBarry Smith   PC_FieldSplitLink ilink = jac->head,next;
845574deadeSBarry Smith 
846574deadeSBarry Smith   PetscFunctionBegin;
847574deadeSBarry Smith   ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr);
848574deadeSBarry Smith   while (ilink) {
8496bf464f9SBarry Smith     ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr);
850574deadeSBarry Smith     next = ilink->next;
851574deadeSBarry Smith     ierr = PetscFree(ilink->splitname);CHKERRQ(ierr);
852574deadeSBarry Smith     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
853574deadeSBarry Smith     ierr = PetscFree(ilink);CHKERRQ(ierr);
854574deadeSBarry Smith     ilink = next;
855574deadeSBarry Smith   }
856574deadeSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
857c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
858d88c663fSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","",PETSC_NULL);CHKERRQ(ierr);
859d88c663fSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","",PETSC_NULL);CHKERRQ(ierr);
860d88c663fSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","",PETSC_NULL);CHKERRQ(ierr);
861d88c663fSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","",PETSC_NULL);CHKERRQ(ierr);
862d88c663fSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","",PETSC_NULL);CHKERRQ(ierr);
8630971522cSBarry Smith   PetscFunctionReturn(0);
8640971522cSBarry Smith }
8650971522cSBarry Smith 
8660971522cSBarry Smith #undef __FUNCT__
8670971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
8680971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
8690971522cSBarry Smith {
8701b9fc7fcSBarry Smith   PetscErrorCode  ierr;
8716c924f48SJed Brown   PetscInt        bs;
872bc59fbc5SBarry Smith   PetscBool       flg,stokes = PETSC_FALSE;
8739dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
8743b224e63SBarry Smith   PCCompositeType ctype;
8751b9fc7fcSBarry Smith 
8760971522cSBarry Smith   PetscFunctionBegin;
8771b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
878acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr);
87951f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
88051f519a2SBarry Smith   if (flg) {
88151f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
88251f519a2SBarry Smith   }
883704ba839SBarry Smith 
884435f959eSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
885c0adfefeSBarry Smith   if (stokes) {
886c0adfefeSBarry Smith     ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
887c0adfefeSBarry Smith     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
888c0adfefeSBarry Smith   }
889c0adfefeSBarry Smith 
8903b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
8913b224e63SBarry Smith   if (flg) {
8923b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
8933b224e63SBarry Smith   }
894d32f9abdSBarry Smith 
895c30613efSMatthew Knepley   /* Only setup fields once */
896c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
897d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
898d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
8996c924f48SJed Brown     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
9006c924f48SJed Brown     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
901d32f9abdSBarry Smith   }
902c5d2311dSJed Brown   if (jac->type == PC_COMPOSITE_SCHUR) {
903c5d2311dSJed Brown     ierr = PetscOptionsEnum("-pc_fieldsplit_schur_factorization_type","Factorization to use","None",PCFieldSplitSchurFactorizationTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,PETSC_NULL);CHKERRQ(ierr);
904084e4875SJed 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);
905c5d2311dSJed Brown   }
9061b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
9070971522cSBarry Smith   PetscFunctionReturn(0);
9080971522cSBarry Smith }
9090971522cSBarry Smith 
9100971522cSBarry Smith /*------------------------------------------------------------------------------------*/
9110971522cSBarry Smith 
9120971522cSBarry Smith EXTERN_C_BEGIN
9130971522cSBarry Smith #undef __FUNCT__
9140971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
9157087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields)
9160971522cSBarry Smith {
91797bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
9180971522cSBarry Smith   PetscErrorCode    ierr;
9195a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
92069a612a9SBarry Smith   char              prefix[128];
92151f519a2SBarry Smith   PetscInt          i;
9220971522cSBarry Smith 
9230971522cSBarry Smith   PetscFunctionBegin;
9246c924f48SJed Brown   if (jac->splitdefined) {
9256c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
9266c924f48SJed Brown     PetscFunctionReturn(0);
9276c924f48SJed Brown   }
92851f519a2SBarry Smith   for (i=0; i<n; i++) {
929e32f2f54SBarry 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);
930e32f2f54SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
93151f519a2SBarry Smith   }
932704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
933a04f6461SBarry Smith   if (splitname) {
934db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
935a04f6461SBarry Smith   } else {
936a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
937a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
938a04f6461SBarry Smith   }
939704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
9405a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
9415a9f2f41SSatish Balay   ilink->nfields = n;
9425a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
9437adad957SLisandro Dalcin   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
9441cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
9455a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
9469005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
94769a612a9SBarry Smith 
948a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
9495a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
9500971522cSBarry Smith 
9510971522cSBarry Smith   if (!next) {
9525a9f2f41SSatish Balay     jac->head       = ilink;
95351f519a2SBarry Smith     ilink->previous = PETSC_NULL;
9540971522cSBarry Smith   } else {
9550971522cSBarry Smith     while (next->next) {
9560971522cSBarry Smith       next = next->next;
9570971522cSBarry Smith     }
9585a9f2f41SSatish Balay     next->next      = ilink;
95951f519a2SBarry Smith     ilink->previous = next;
9600971522cSBarry Smith   }
9610971522cSBarry Smith   jac->nsplits++;
9620971522cSBarry Smith   PetscFunctionReturn(0);
9630971522cSBarry Smith }
9640971522cSBarry Smith EXTERN_C_END
9650971522cSBarry Smith 
966e69d4d44SBarry Smith EXTERN_C_BEGIN
967e69d4d44SBarry Smith #undef __FUNCT__
968e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
9697087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
970e69d4d44SBarry Smith {
971e69d4d44SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
972e69d4d44SBarry Smith   PetscErrorCode ierr;
973e69d4d44SBarry Smith 
974e69d4d44SBarry Smith   PetscFunctionBegin;
975519d70e2SJed Brown   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
976e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
977e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
97813e0d083SBarry Smith   if (n) *n = jac->nsplits;
979e69d4d44SBarry Smith   PetscFunctionReturn(0);
980e69d4d44SBarry Smith }
981e69d4d44SBarry Smith EXTERN_C_END
9820971522cSBarry Smith 
9830971522cSBarry Smith EXTERN_C_BEGIN
9840971522cSBarry Smith #undef __FUNCT__
98569a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
9867087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
9870971522cSBarry Smith {
9880971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
9890971522cSBarry Smith   PetscErrorCode    ierr;
9900971522cSBarry Smith   PetscInt          cnt = 0;
9915a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
9920971522cSBarry Smith 
9930971522cSBarry Smith   PetscFunctionBegin;
9945d480477SMatthew G Knepley   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
9955a9f2f41SSatish Balay   while (ilink) {
9965a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
9975a9f2f41SSatish Balay     ilink = ilink->next;
9980971522cSBarry Smith   }
9995d480477SMatthew 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);
100013e0d083SBarry Smith   if (n) *n = jac->nsplits;
10010971522cSBarry Smith   PetscFunctionReturn(0);
10020971522cSBarry Smith }
10030971522cSBarry Smith EXTERN_C_END
10040971522cSBarry Smith 
1005704ba839SBarry Smith EXTERN_C_BEGIN
1006704ba839SBarry Smith #undef __FUNCT__
1007704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
10087087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
1009704ba839SBarry Smith {
1010704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1011704ba839SBarry Smith   PetscErrorCode    ierr;
1012704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
1013704ba839SBarry Smith   char              prefix[128];
1014704ba839SBarry Smith 
1015704ba839SBarry Smith   PetscFunctionBegin;
10166c924f48SJed Brown   if (jac->splitdefined) {
10176c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
10186c924f48SJed Brown     PetscFunctionReturn(0);
10196c924f48SJed Brown   }
102016913363SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
1021a04f6461SBarry Smith   if (splitname) {
1022db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1023a04f6461SBarry Smith   } else {
1024a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
1025d88c663fSJed Brown     ierr = PetscSNPrintf(ilink->splitname,3,"%D",jac->nsplits);CHKERRQ(ierr);
1026a04f6461SBarry Smith   }
10271ab39975SBarry Smith   ilink->is      = is;
10281ab39975SBarry Smith   ierr           = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1029704ba839SBarry Smith   ilink->next    = PETSC_NULL;
1030704ba839SBarry Smith   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
10311cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1032704ba839SBarry Smith   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
10339005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
1034704ba839SBarry Smith 
1035a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
1036704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
1037704ba839SBarry Smith 
1038704ba839SBarry Smith   if (!next) {
1039704ba839SBarry Smith     jac->head       = ilink;
1040704ba839SBarry Smith     ilink->previous = PETSC_NULL;
1041704ba839SBarry Smith   } else {
1042704ba839SBarry Smith     while (next->next) {
1043704ba839SBarry Smith       next = next->next;
1044704ba839SBarry Smith     }
1045704ba839SBarry Smith     next->next      = ilink;
1046704ba839SBarry Smith     ilink->previous = next;
1047704ba839SBarry Smith   }
1048704ba839SBarry Smith   jac->nsplits++;
1049704ba839SBarry Smith 
1050704ba839SBarry Smith   PetscFunctionReturn(0);
1051704ba839SBarry Smith }
1052704ba839SBarry Smith EXTERN_C_END
1053704ba839SBarry Smith 
10540971522cSBarry Smith #undef __FUNCT__
10550971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
10560971522cSBarry Smith /*@
10570971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
10580971522cSBarry Smith 
1059ad4df100SBarry Smith     Logically Collective on PC
10600971522cSBarry Smith 
10610971522cSBarry Smith     Input Parameters:
10620971522cSBarry Smith +   pc  - the preconditioner context
1063a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
10640971522cSBarry Smith .   n - the number of fields in this split
1065db4c96c1SJed Brown -   fields - the fields in this split
10660971522cSBarry Smith 
10670971522cSBarry Smith     Level: intermediate
10680971522cSBarry Smith 
1069d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
1070d32f9abdSBarry Smith 
1071d32f9abdSBarry Smith      The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block
1072d32f9abdSBarry 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
1073d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1074d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
1075d32f9abdSBarry Smith 
1076db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
1077db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1078db4c96c1SJed Brown 
1079d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
10800971522cSBarry Smith 
10810971522cSBarry Smith @*/
10827087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields)
10830971522cSBarry Smith {
10844ac538c5SBarry Smith   PetscErrorCode ierr;
10850971522cSBarry Smith 
10860971522cSBarry Smith   PetscFunctionBegin;
10870700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1088db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1089db4c96c1SJed Brown   if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1090db4c96c1SJed Brown   PetscValidIntPointer(fields,3);
10914ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *),(pc,splitname,n,fields));CHKERRQ(ierr);
10920971522cSBarry Smith   PetscFunctionReturn(0);
10930971522cSBarry Smith }
10940971522cSBarry Smith 
10950971522cSBarry Smith #undef __FUNCT__
1096704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
1097704ba839SBarry Smith /*@
1098704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
1099704ba839SBarry Smith 
1100ad4df100SBarry Smith     Logically Collective on PC
1101704ba839SBarry Smith 
1102704ba839SBarry Smith     Input Parameters:
1103704ba839SBarry Smith +   pc  - the preconditioner context
1104a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
1105db4c96c1SJed Brown -   is - the index set that defines the vector elements in this field
1106704ba839SBarry Smith 
1107d32f9abdSBarry Smith 
1108a6ffb8dbSJed Brown     Notes:
1109a6ffb8dbSJed Brown     Use PCFieldSplitSetFields(), for fields defined by strided types.
1110a6ffb8dbSJed Brown 
1111db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
1112db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1113d32f9abdSBarry Smith 
1114704ba839SBarry Smith     Level: intermediate
1115704ba839SBarry Smith 
1116704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
1117704ba839SBarry Smith 
1118704ba839SBarry Smith @*/
11197087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1120704ba839SBarry Smith {
11214ac538c5SBarry Smith   PetscErrorCode ierr;
1122704ba839SBarry Smith 
1123704ba839SBarry Smith   PetscFunctionBegin;
11240700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1125d88c663fSJed Brown   if (splitname) PetscValidCharPointer(splitname,2);
1126db4c96c1SJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,3);
11274ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr);
1128704ba839SBarry Smith   PetscFunctionReturn(0);
1129704ba839SBarry Smith }
1130704ba839SBarry Smith 
1131704ba839SBarry Smith #undef __FUNCT__
113257a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS"
113357a9adfeSMatthew G Knepley /*@
113457a9adfeSMatthew G Knepley     PCFieldSplitGetIS - Retrieves the elements for a field as an IS
113557a9adfeSMatthew G Knepley 
113657a9adfeSMatthew G Knepley     Logically Collective on PC
113757a9adfeSMatthew G Knepley 
113857a9adfeSMatthew G Knepley     Input Parameters:
113957a9adfeSMatthew G Knepley +   pc  - the preconditioner context
114057a9adfeSMatthew G Knepley -   splitname - name of this split
114157a9adfeSMatthew G Knepley 
114257a9adfeSMatthew G Knepley     Output Parameter:
114357a9adfeSMatthew G Knepley -   is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found
114457a9adfeSMatthew G Knepley 
114557a9adfeSMatthew G Knepley     Level: intermediate
114657a9adfeSMatthew G Knepley 
114757a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS()
114857a9adfeSMatthew G Knepley 
114957a9adfeSMatthew G Knepley @*/
115057a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
115157a9adfeSMatthew G Knepley {
115257a9adfeSMatthew G Knepley   PetscErrorCode ierr;
115357a9adfeSMatthew G Knepley 
115457a9adfeSMatthew G Knepley   PetscFunctionBegin;
115557a9adfeSMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
115657a9adfeSMatthew G Knepley   PetscValidCharPointer(splitname,2);
115757a9adfeSMatthew G Knepley   PetscValidPointer(is,3);
115857a9adfeSMatthew G Knepley   {
115957a9adfeSMatthew G Knepley     PC_FieldSplit    *jac   = (PC_FieldSplit *) pc->data;
116057a9adfeSMatthew G Knepley     PC_FieldSplitLink ilink = jac->head;
116157a9adfeSMatthew G Knepley     PetscBool         found;
116257a9adfeSMatthew G Knepley 
116357a9adfeSMatthew G Knepley     *is = PETSC_NULL;
116457a9adfeSMatthew G Knepley     while(ilink) {
116557a9adfeSMatthew G Knepley       ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr);
116657a9adfeSMatthew G Knepley       if (found) {
116757a9adfeSMatthew G Knepley         *is = ilink->is;
116857a9adfeSMatthew G Knepley         break;
116957a9adfeSMatthew G Knepley       }
117057a9adfeSMatthew G Knepley       ilink = ilink->next;
117157a9adfeSMatthew G Knepley     }
117257a9adfeSMatthew G Knepley   }
117357a9adfeSMatthew G Knepley   PetscFunctionReturn(0);
117457a9adfeSMatthew G Knepley }
117557a9adfeSMatthew G Knepley 
117657a9adfeSMatthew G Knepley #undef __FUNCT__
117751f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
117851f519a2SBarry Smith /*@
117951f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
118051f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
118151f519a2SBarry Smith 
1182ad4df100SBarry Smith     Logically Collective on PC
118351f519a2SBarry Smith 
118451f519a2SBarry Smith     Input Parameters:
118551f519a2SBarry Smith +   pc  - the preconditioner context
118651f519a2SBarry Smith -   bs - the block size
118751f519a2SBarry Smith 
118851f519a2SBarry Smith     Level: intermediate
118951f519a2SBarry Smith 
119051f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
119151f519a2SBarry Smith 
119251f519a2SBarry Smith @*/
11937087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
119451f519a2SBarry Smith {
11954ac538c5SBarry Smith   PetscErrorCode ierr;
119651f519a2SBarry Smith 
119751f519a2SBarry Smith   PetscFunctionBegin;
11980700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1199c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,bs,2);
12004ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr);
120151f519a2SBarry Smith   PetscFunctionReturn(0);
120251f519a2SBarry Smith }
120351f519a2SBarry Smith 
120451f519a2SBarry Smith #undef __FUNCT__
120569a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
12060971522cSBarry Smith /*@C
120769a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
12080971522cSBarry Smith 
120969a612a9SBarry Smith    Collective on KSP
12100971522cSBarry Smith 
12110971522cSBarry Smith    Input Parameter:
12120971522cSBarry Smith .  pc - the preconditioner context
12130971522cSBarry Smith 
12140971522cSBarry Smith    Output Parameters:
121513e0d083SBarry Smith +  n - the number of splits
121669a612a9SBarry Smith -  pc - the array of KSP contexts
12170971522cSBarry Smith 
12180971522cSBarry Smith    Note:
1219d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
1220d32f9abdSBarry Smith    (not the KSP just the array that contains them).
12210971522cSBarry Smith 
122269a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
12230971522cSBarry Smith 
12240971522cSBarry Smith    Level: advanced
12250971522cSBarry Smith 
12260971522cSBarry Smith .seealso: PCFIELDSPLIT
12270971522cSBarry Smith @*/
12287087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
12290971522cSBarry Smith {
12304ac538c5SBarry Smith   PetscErrorCode ierr;
12310971522cSBarry Smith 
12320971522cSBarry Smith   PetscFunctionBegin;
12330700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
123413e0d083SBarry Smith   if (n) PetscValidIntPointer(n,2);
12354ac538c5SBarry Smith   ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
12360971522cSBarry Smith   PetscFunctionReturn(0);
12370971522cSBarry Smith }
12380971522cSBarry Smith 
1239e69d4d44SBarry Smith #undef __FUNCT__
1240e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition"
1241e69d4d44SBarry Smith /*@
1242e69d4d44SBarry Smith     PCFieldSplitSchurPrecondition -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
1243a04f6461SBarry Smith       A11 matrix. Otherwise no preconditioner is used.
1244e69d4d44SBarry Smith 
1245e69d4d44SBarry Smith     Collective on PC
1246e69d4d44SBarry Smith 
1247e69d4d44SBarry Smith     Input Parameters:
1248e69d4d44SBarry Smith +   pc  - the preconditioner context
1249fd1303e9SJungho Lee .   ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default
1250084e4875SJed Brown -   userpre - matrix to use for preconditioning, or PETSC_NULL
1251084e4875SJed Brown 
1252e69d4d44SBarry Smith     Options Database:
1253084e4875SJed Brown .     -pc_fieldsplit_schur_precondition <self,user,diag> default is diag
1254e69d4d44SBarry Smith 
1255fd1303e9SJungho Lee     Notes:
1256fd1303e9SJungho Lee $    If ptype is
1257fd1303e9SJungho Lee $        user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument
1258fd1303e9SJungho Lee $             to this function).
1259fd1303e9SJungho Lee $        diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original
1260fd1303e9SJungho Lee $             matrix associated with the Schur complement (i.e. A11)
1261fd1303e9SJungho Lee $        self the preconditioner for the Schur complement is generated from the Schur complement matrix itself:
1262fd1303e9SJungho Lee $             The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC
1263fd1303e9SJungho Lee $             preconditioner
1264fd1303e9SJungho Lee 
1265fd1303e9SJungho Lee      When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense
1266fd1303e9SJungho Lee     with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
1267fd1303e9SJungho Lee     -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement.
1268fd1303e9SJungho Lee 
1269fd1303e9SJungho Lee     Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in
1270fd1303e9SJungho Lee      the name since it sets a proceedure to use.
1271fd1303e9SJungho Lee 
1272e69d4d44SBarry Smith     Level: intermediate
1273e69d4d44SBarry Smith 
1274fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC
1275e69d4d44SBarry Smith 
1276e69d4d44SBarry Smith @*/
12777087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1278e69d4d44SBarry Smith {
12794ac538c5SBarry Smith   PetscErrorCode ierr;
1280e69d4d44SBarry Smith 
1281e69d4d44SBarry Smith   PetscFunctionBegin;
12820700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12834ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr);
1284e69d4d44SBarry Smith   PetscFunctionReturn(0);
1285e69d4d44SBarry Smith }
1286e69d4d44SBarry Smith 
1287e69d4d44SBarry Smith EXTERN_C_BEGIN
1288e69d4d44SBarry Smith #undef __FUNCT__
1289e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
12907087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1291e69d4d44SBarry Smith {
1292e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1293084e4875SJed Brown   PetscErrorCode  ierr;
1294e69d4d44SBarry Smith 
1295e69d4d44SBarry Smith   PetscFunctionBegin;
1296084e4875SJed Brown   jac->schurpre = ptype;
1297084e4875SJed Brown   if (pre) {
12986bf464f9SBarry Smith     ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
1299084e4875SJed Brown     jac->schur_user = pre;
1300084e4875SJed Brown     ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1301084e4875SJed Brown   }
1302e69d4d44SBarry Smith   PetscFunctionReturn(0);
1303e69d4d44SBarry Smith }
1304e69d4d44SBarry Smith EXTERN_C_END
1305e69d4d44SBarry Smith 
130630ad9308SMatthew Knepley #undef __FUNCT__
130730ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
130830ad9308SMatthew Knepley /*@C
13098c03b21aSDmitry Karpeev    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
131030ad9308SMatthew Knepley 
131130ad9308SMatthew Knepley    Collective on KSP
131230ad9308SMatthew Knepley 
131330ad9308SMatthew Knepley    Input Parameter:
131430ad9308SMatthew Knepley .  pc - the preconditioner context
131530ad9308SMatthew Knepley 
131630ad9308SMatthew Knepley    Output Parameters:
1317a04f6461SBarry Smith +  A00 - the (0,0) block
1318a04f6461SBarry Smith .  A01 - the (0,1) block
1319a04f6461SBarry Smith .  A10 - the (1,0) block
1320a04f6461SBarry Smith -  A11 - the (1,1) block
132130ad9308SMatthew Knepley 
132230ad9308SMatthew Knepley    Level: advanced
132330ad9308SMatthew Knepley 
132430ad9308SMatthew Knepley .seealso: PCFIELDSPLIT
132530ad9308SMatthew Knepley @*/
1326a04f6461SBarry Smith PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
132730ad9308SMatthew Knepley {
132830ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit *) pc->data;
132930ad9308SMatthew Knepley 
133030ad9308SMatthew Knepley   PetscFunctionBegin;
13310700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1332c1235816SBarry Smith   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
1333a04f6461SBarry Smith   if (A00) *A00 = jac->pmat[0];
1334a04f6461SBarry Smith   if (A01) *A01 = jac->B;
1335a04f6461SBarry Smith   if (A10) *A10 = jac->C;
1336a04f6461SBarry Smith   if (A11) *A11 = jac->pmat[1];
133730ad9308SMatthew Knepley   PetscFunctionReturn(0);
133830ad9308SMatthew Knepley }
133930ad9308SMatthew Knepley 
134079416396SBarry Smith EXTERN_C_BEGIN
134179416396SBarry Smith #undef __FUNCT__
134279416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
13437087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
134479416396SBarry Smith {
134579416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1346e69d4d44SBarry Smith   PetscErrorCode ierr;
134779416396SBarry Smith 
134879416396SBarry Smith   PetscFunctionBegin;
134979416396SBarry Smith   jac->type = type;
13503b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
13513b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
13523b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
1353e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1354e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1355e69d4d44SBarry Smith 
13563b224e63SBarry Smith   } else {
13573b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
13583b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
1359e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
13609e7d6b0aSBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr);
13613b224e63SBarry Smith   }
136279416396SBarry Smith   PetscFunctionReturn(0);
136379416396SBarry Smith }
136479416396SBarry Smith EXTERN_C_END
136579416396SBarry Smith 
136651f519a2SBarry Smith EXTERN_C_BEGIN
136751f519a2SBarry Smith #undef __FUNCT__
136851f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
13697087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
137051f519a2SBarry Smith {
137151f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
137251f519a2SBarry Smith 
137351f519a2SBarry Smith   PetscFunctionBegin;
137465e19b50SBarry Smith   if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
137565e19b50SBarry 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);
137651f519a2SBarry Smith   jac->bs = bs;
137751f519a2SBarry Smith   PetscFunctionReturn(0);
137851f519a2SBarry Smith }
137951f519a2SBarry Smith EXTERN_C_END
138051f519a2SBarry Smith 
138179416396SBarry Smith #undef __FUNCT__
138279416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
1383bc08b0f1SBarry Smith /*@
138479416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
138579416396SBarry Smith 
138679416396SBarry Smith    Collective on PC
138779416396SBarry Smith 
138879416396SBarry Smith    Input Parameter:
138979416396SBarry Smith .  pc - the preconditioner context
139081540f2fSBarry Smith .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
139179416396SBarry Smith 
139279416396SBarry Smith    Options Database Key:
1393a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
139479416396SBarry Smith 
139579416396SBarry Smith    Level: Developer
139679416396SBarry Smith 
139779416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
139879416396SBarry Smith 
139979416396SBarry Smith .seealso: PCCompositeSetType()
140079416396SBarry Smith 
140179416396SBarry Smith @*/
14027087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
140379416396SBarry Smith {
14044ac538c5SBarry Smith   PetscErrorCode ierr;
140579416396SBarry Smith 
140679416396SBarry Smith   PetscFunctionBegin;
14070700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
14084ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr);
140979416396SBarry Smith   PetscFunctionReturn(0);
141079416396SBarry Smith }
141179416396SBarry Smith 
14120971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
14130971522cSBarry Smith /*MC
1414a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
1415a04f6461SBarry Smith                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.
14160971522cSBarry Smith 
1417edf189efSBarry Smith      To set options on the solvers for each block append -fieldsplit_ to all the PC
1418edf189efSBarry Smith         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
14190971522cSBarry Smith 
1420a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
142169a612a9SBarry Smith          and set the options directly on the resulting KSP object
14220971522cSBarry Smith 
14230971522cSBarry Smith    Level: intermediate
14240971522cSBarry Smith 
142579416396SBarry Smith    Options Database Keys:
142681540f2fSBarry Smith +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
142781540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
142881540f2fSBarry Smith                               been supplied explicitly by -pc_fieldsplit_%d_fields
142981540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
143081540f2fSBarry Smith .   -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative>
1431e69d4d44SBarry Smith .   -pc_fieldsplit_schur_precondition <true,false> default is true
1432435f959eSBarry Smith .   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver
143343890ad2SJed Brown .   -fieldsplit_NAME_ksp_* - control inner linear solver, NAME is a sequential integer if unspecified, otherwise use name provided in PCFieldSplitSetIS() or the name associated with the field in the DM passed to PCSetDM() (or a higher level function)
143443890ad2SJed Brown -   -fieldsplit_NAME_pc_* - control inner preconditioner, NAME has same semantics as above
143579416396SBarry Smith 
143643890ad2SJed Brown    Notes:
143743890ad2SJed Brown     Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1438d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
1439d32f9abdSBarry Smith 
1440d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
1441d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1442d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
1443d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
1444d32f9abdSBarry Smith 
1445a04f6461SBarry Smith      For the Schur complement preconditioner if J = ( A00 A01 )
1446a04f6461SBarry Smith                                                     ( A10 A11 )
1447e69d4d44SBarry Smith      the preconditioner is
1448a04f6461SBarry Smith               (I   -B ksp(A00)) ( inv(A00)   0  (I         0  )
1449a04f6461SBarry Smith               (0    I       ) (   0    ksp(S) ) (-A10 ksp(A00) I  )
1450a04f6461SBarry Smith      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of
1451a04f6461SBarry 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).
1452a04f6461SBarry Smith      For PCFieldSplitGetKSP() when field number is
1453edf189efSBarry Smith      0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
1454a04f6461SBarry Smith      A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
14557e8cb189SBarry Smith      option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc
1456e69d4d44SBarry Smith 
1457edf189efSBarry Smith      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
1458edf189efSBarry Smith      is used automatically for a second block.
1459edf189efSBarry Smith 
1460ff218e97SBarry Smith      The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1.
1461ff218e97SBarry Smith      Generally it should be used with the AIJ format.
1462ff218e97SBarry Smith 
1463ff218e97SBarry Smith      The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see,
1464ff218e97SBarry Smith      for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT
1465ff218e97SBarry Smith      inside a smoother resulting in "Distributive Smoothers".
14660716a85fSBarry Smith 
1467a541d17aSBarry Smith    Concepts: physics based preconditioners, block preconditioners
14680971522cSBarry Smith 
14697e8cb189SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
1470e69d4d44SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
14710971522cSBarry Smith M*/
14720971522cSBarry Smith 
14730971522cSBarry Smith EXTERN_C_BEGIN
14740971522cSBarry Smith #undef __FUNCT__
14750971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
14767087cfbeSBarry Smith PetscErrorCode  PCCreate_FieldSplit(PC pc)
14770971522cSBarry Smith {
14780971522cSBarry Smith   PetscErrorCode ierr;
14790971522cSBarry Smith   PC_FieldSplit  *jac;
14800971522cSBarry Smith 
14810971522cSBarry Smith   PetscFunctionBegin;
148238f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
14830971522cSBarry Smith   jac->bs        = -1;
14840971522cSBarry Smith   jac->nsplits   = 0;
14853e197d65SBarry Smith   jac->type      = PC_COMPOSITE_MULTIPLICATIVE;
1486e6cab6aaSJed Brown   jac->schurpre  = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
1487c5d2311dSJed Brown   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL;
148851f519a2SBarry Smith 
14890971522cSBarry Smith   pc->data     = (void*)jac;
14900971522cSBarry Smith 
14910971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
1492421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
14930971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
1494574deadeSBarry Smith   pc->ops->reset             = PCReset_FieldSplit;
14950971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
14960971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
14970971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
14980971522cSBarry Smith   pc->ops->applyrichardson   = 0;
14990971522cSBarry Smith 
150069a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
150169a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
15020971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
15030971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1504704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
1505704ba839SBarry Smith                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
150679416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
150779416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
150851f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
150951f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
15100971522cSBarry Smith   PetscFunctionReturn(0);
15110971522cSBarry Smith }
15120971522cSBarry Smith EXTERN_C_END
15130971522cSBarry Smith 
15140971522cSBarry Smith 
1515a541d17aSBarry Smith 
1516