xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision d282228771c9a5fecbd3a3cc324bb1794e538e5d)
1dba47a55SKris Buschelman 
2c6db04a5SJed Brown #include <private/pcimpl.h>     /*I "petscpc.h" I*/
33c48a1e8SJed Brown #include <petscdmcomposite.h>   /*I "petscdmcomposite.h" I*/
4*d2822287SMatthew G Knepley #include <petscdmcomplex.h>
50971522cSBarry Smith 
6f7c28744SJed Brown const char *const PCFieldSplitSchurPreTypes[] = {"SELF","DIAG","USER","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0};
7f7c28744SJed Brown const char *const PCFieldSplitSchurFactorizationTypes[] = {"DIAG","LOWER","UPPER","FULL","PCFieldSplitSchurFactorizationType","PC_FIELDSPLIT_SCHUR_FACTORIZATION_",0};
8c5d2311dSJed Brown 
9c5d2311dSJed Brown typedef enum {
10c5d2311dSJed Brown   PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG,
11c5d2311dSJed Brown   PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER,
12c5d2311dSJed Brown   PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER,
13c5d2311dSJed Brown   PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL
14c5d2311dSJed Brown } PCFieldSplitSchurFactorizationType;
15084e4875SJed Brown 
160971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
170971522cSBarry Smith struct _PC_FieldSplitLink {
1869a612a9SBarry Smith   KSP               ksp;
190971522cSBarry Smith   Vec               x,y;
20db4c96c1SJed Brown   char              *splitname;
210971522cSBarry Smith   PetscInt          nfields;
220971522cSBarry Smith   PetscInt          *fields;
231b9fc7fcSBarry Smith   VecScatter        sctx;
244aa3045dSJed Brown   IS                is;
2551f519a2SBarry Smith   PC_FieldSplitLink next,previous;
260971522cSBarry Smith };
270971522cSBarry Smith 
280971522cSBarry Smith typedef struct {
2968dd23aaSBarry Smith   PCCompositeType                    type;
30ace3abfcSBarry Smith   PetscBool                          defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */
31ace3abfcSBarry Smith   PetscBool                          splitdefined; /* Flag is set after the splits have been defined, to prevent more splits from being added */
32ace3abfcSBarry Smith   PetscBool                          realdiagonal; /* Flag to use the diagonal blocks of mat preconditioned by pmat, instead of just pmat */
3330ad9308SMatthew Knepley   PetscInt                           bs;           /* Block size for IS and Mat structures */
3430ad9308SMatthew Knepley   PetscInt                           nsplits;      /* Number of field divisions defined */
3579416396SBarry Smith   Vec                                *x,*y,w1,w2;
36519d70e2SJed Brown   Mat                                *mat;         /* The diagonal block for each split */
37519d70e2SJed Brown   Mat                                *pmat;        /* The preconditioning diagonal block for each split */
3830ad9308SMatthew Knepley   Mat                                *Afield;      /* The rows of the matrix associated with each split */
39ace3abfcSBarry Smith   PetscBool                          issetup;
4030ad9308SMatthew Knepley   /* Only used when Schur complement preconditioning is used */
4130ad9308SMatthew Knepley   Mat                                B;            /* The (0,1) block */
4230ad9308SMatthew Knepley   Mat                                C;            /* The (1,0) block */
43a04f6461SBarry Smith   Mat                                schur;        /* The Schur complement S = A11 - A10 A00^{-1} A01 */
44084e4875SJed Brown   Mat                                schur_user;   /* User-provided preconditioning matrix for the Schur complement */
45084e4875SJed Brown   PCFieldSplitSchurPreType           schurpre; /* Determines which preconditioning matrix is used for the Schur complement */
46c5d2311dSJed Brown   PCFieldSplitSchurFactorizationType schurfactorization;
4730ad9308SMatthew Knepley   KSP                                kspschur;     /* The solver for S */
4897bbdb24SBarry Smith   PC_FieldSplitLink                  head;
4963ec74ffSBarry Smith   PetscBool                          reset;         /* indicates PCReset() has been last called on this object, hack */
50c1570756SJed Brown   PetscBool                          suboptionsset; /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */
510971522cSBarry Smith } PC_FieldSplit;
520971522cSBarry Smith 
5316913363SBarry Smith /*
5416913363SBarry Smith     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
5516913363SBarry Smith    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
5616913363SBarry Smith    PC you could change this.
5716913363SBarry Smith */
58084e4875SJed Brown 
59e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it.  This way the
60084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */
61084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac)
62084e4875SJed Brown {
63084e4875SJed Brown   switch (jac->schurpre) {
64084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur;
65084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1];
66084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */
67084e4875SJed Brown     default:
68084e4875SJed Brown       return jac->schur_user ? jac->schur_user : jac->pmat[1];
69084e4875SJed Brown   }
70084e4875SJed Brown }
71084e4875SJed Brown 
72084e4875SJed Brown 
730971522cSBarry Smith #undef __FUNCT__
740971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
750971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
760971522cSBarry Smith {
770971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
780971522cSBarry Smith   PetscErrorCode    ierr;
79ace3abfcSBarry Smith   PetscBool         iascii;
800971522cSBarry Smith   PetscInt          i,j;
815a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
820971522cSBarry Smith 
830971522cSBarry Smith   PetscFunctionBegin;
842692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
850971522cSBarry Smith   if (iascii) {
8651f519a2SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
8769a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
880971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
890971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
901ab39975SBarry Smith       if (ilink->fields) {
910971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
9279416396SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
935a9f2f41SSatish Balay 	for (j=0; j<ilink->nfields; j++) {
9479416396SBarry Smith 	  if (j > 0) {
9579416396SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
9679416396SBarry Smith 	  }
975a9f2f41SSatish Balay 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
980971522cSBarry Smith 	}
990971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
10079416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1011ab39975SBarry Smith       } else {
1021ab39975SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1031ab39975SBarry Smith       }
1045a9f2f41SSatish Balay       ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
1055a9f2f41SSatish Balay       ilink = ilink->next;
1060971522cSBarry Smith     }
1070971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1080971522cSBarry Smith   } else {
10965e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1100971522cSBarry Smith   }
1110971522cSBarry Smith   PetscFunctionReturn(0);
1120971522cSBarry Smith }
1130971522cSBarry Smith 
1140971522cSBarry Smith #undef __FUNCT__
1153b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur"
1163b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
1173b224e63SBarry Smith {
1183b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1193b224e63SBarry Smith   PetscErrorCode    ierr;
120ace3abfcSBarry Smith   PetscBool         iascii;
1213b224e63SBarry Smith   PetscInt          i,j;
1223b224e63SBarry Smith   PC_FieldSplitLink ilink = jac->head;
1233b224e63SBarry Smith   KSP               ksp;
1243b224e63SBarry Smith 
1253b224e63SBarry Smith   PetscFunctionBegin;
1262692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1273b224e63SBarry Smith   if (iascii) {
128c5d2311dSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactorizationTypes[jac->schurfactorization]);CHKERRQ(ierr);
1293b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Split info:\n");CHKERRQ(ierr);
1303b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1313b224e63SBarry Smith     for (i=0; i<jac->nsplits; i++) {
1323b224e63SBarry Smith       if (ilink->fields) {
1333b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
1343b224e63SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
1353b224e63SBarry Smith 	for (j=0; j<ilink->nfields; j++) {
1363b224e63SBarry Smith 	  if (j > 0) {
1373b224e63SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
1383b224e63SBarry Smith 	  }
1393b224e63SBarry Smith 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
1403b224e63SBarry Smith 	}
1413b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1423b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1433b224e63SBarry Smith       } else {
1443b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1453b224e63SBarry Smith       }
1463b224e63SBarry Smith       ilink = ilink->next;
1473b224e63SBarry Smith     }
148435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block \n");CHKERRQ(ierr);
1493b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15012cae6f2SJed Brown     if (jac->schur) {
15112cae6f2SJed Brown       ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
1523b224e63SBarry Smith       ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
15312cae6f2SJed Brown     } else {
15412cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
15512cae6f2SJed Brown     }
1563b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
157435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr);
1583b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15912cae6f2SJed Brown     if (jac->kspschur) {
1603b224e63SBarry Smith       ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
16112cae6f2SJed Brown     } else {
16212cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
16312cae6f2SJed Brown     }
1643b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1653b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1663b224e63SBarry Smith   } else {
16765e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1683b224e63SBarry Smith   }
1693b224e63SBarry Smith   PetscFunctionReturn(0);
1703b224e63SBarry Smith }
1713b224e63SBarry Smith 
1723b224e63SBarry Smith #undef __FUNCT__
1736c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private"
1746c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */
1756c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc)
1766c924f48SJed Brown {
1776c924f48SJed Brown   PetscErrorCode ierr;
1786c924f48SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1796c924f48SJed Brown   PetscInt       i,nfields,*ifields;
180ace3abfcSBarry Smith   PetscBool      flg;
1816c924f48SJed Brown   char           optionname[128],splitname[8];
1826c924f48SJed Brown 
1836c924f48SJed Brown   PetscFunctionBegin;
1846c924f48SJed Brown   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr);
1856c924f48SJed Brown   for (i=0,flg=PETSC_TRUE; ; i++) {
1866c924f48SJed Brown     ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
1876c924f48SJed Brown     ierr = PetscSNPrintf(optionname,sizeof optionname,"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr);
1886c924f48SJed Brown     nfields = jac->bs;
1896c924f48SJed Brown     ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr);
1906c924f48SJed Brown     if (!flg) break;
191ea6813d4SMatthew G Knepley     if (!nfields) SETERRQ(((PetscObject) pc)->comm,PETSC_ERR_USER,"Cannot list zero fields");
1926c924f48SJed Brown     ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields);CHKERRQ(ierr);
1936c924f48SJed Brown   }
1946c924f48SJed Brown   if (i > 0) {
1956c924f48SJed Brown     /* Makes command-line setting of splits take precedence over setting them in code.
1966c924f48SJed Brown        Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would
1976c924f48SJed Brown        create new splits, which would probably not be what the user wanted. */
1986c924f48SJed Brown     jac->splitdefined = PETSC_TRUE;
1996c924f48SJed Brown   }
2006c924f48SJed Brown   ierr = PetscFree(ifields);CHKERRQ(ierr);
2016c924f48SJed Brown   PetscFunctionReturn(0);
2026c924f48SJed Brown }
2036c924f48SJed Brown 
2046c924f48SJed Brown #undef __FUNCT__
20569a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
20669a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
2070971522cSBarry Smith {
2080971522cSBarry Smith   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
2090971522cSBarry Smith   PetscErrorCode    ierr;
2105a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
2116ce1633cSBarry Smith   PetscBool         flg = PETSC_FALSE,stokes = PETSC_FALSE;
2126c924f48SJed Brown   PetscInt          i;
2130971522cSBarry Smith 
2140971522cSBarry Smith   PetscFunctionBegin;
215d32f9abdSBarry Smith   if (!ilink) {
216d0af7cd3SBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
217d0af7cd3SBarry Smith     if (pc->dm && !stokes) {
218ea6813d4SMatthew G Knepley       PetscBool dmcomposite, dmcomplex;
2198b8307b2SJed Brown       ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr);
220ea6813d4SMatthew G Knepley       ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPLEX,&dmcomplex);CHKERRQ(ierr);
2218b8307b2SJed Brown       if (dmcomposite) {
2228b8307b2SJed Brown         PetscInt nDM;
2238b8307b2SJed Brown         IS       *fields;
2242fa5ba8aSJed Brown         DM       *dms;
2258b8307b2SJed Brown         ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
2268b8307b2SJed Brown         ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr);
2278b8307b2SJed Brown         ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr);
2282fa5ba8aSJed Brown         ierr = PetscMalloc(nDM*sizeof(DM),&dms);CHKERRQ(ierr);
2292fa5ba8aSJed Brown         ierr = DMCompositeGetEntriesArray(pc->dm,dms);CHKERRQ(ierr);
2308b8307b2SJed Brown         for (i=0; i<nDM; i++) {
2312fa5ba8aSJed Brown           char buf[256];
2322fa5ba8aSJed Brown           const char *splitname;
2332fa5ba8aSJed Brown           /* Split naming precedence: object name, prefix, number */
2342fa5ba8aSJed Brown           splitname = ((PetscObject)pc->dm)->name;
2352fa5ba8aSJed Brown           if (!splitname) {
2362fa5ba8aSJed Brown             ierr = PetscObjectGetOptionsPrefix((PetscObject)dms[i],&splitname);CHKERRQ(ierr);
2372fa5ba8aSJed Brown             if (splitname) {
2382fa5ba8aSJed Brown               size_t len;
2392fa5ba8aSJed Brown               ierr = PetscStrncpy(buf,splitname,sizeof buf);CHKERRQ(ierr);
2402fa5ba8aSJed Brown               buf[sizeof buf - 1] = 0;
2412fa5ba8aSJed Brown               ierr = PetscStrlen(buf,&len);CHKERRQ(ierr);
2422fa5ba8aSJed Brown               if (buf[len-1] == '_') buf[len-1] = 0; /* Remove trailing underscore if it was used */
2432fa5ba8aSJed Brown               splitname = buf;
2442fa5ba8aSJed Brown             }
2452fa5ba8aSJed Brown           }
2462fa5ba8aSJed Brown           if (!splitname) {
2472fa5ba8aSJed Brown             ierr = PetscSNPrintf(buf,sizeof buf,"%D",i);CHKERRQ(ierr);
2482fa5ba8aSJed Brown             splitname = buf;
2492fa5ba8aSJed Brown           }
2508b8307b2SJed Brown           ierr = PCFieldSplitSetIS(pc,splitname,fields[i]);CHKERRQ(ierr);
251fcfd50ebSBarry Smith           ierr = ISDestroy(&fields[i]);CHKERRQ(ierr);
2528b8307b2SJed Brown         }
2538b8307b2SJed Brown         ierr = PetscFree(fields);CHKERRQ(ierr);
2542fa5ba8aSJed Brown         for (ilink=jac->head,i=0; ilink; ilink=ilink->next,i++) {
2552fa5ba8aSJed Brown           ierr = KSPSetDM(ilink->ksp,dms[i]);CHKERRQ(ierr);
2562fa5ba8aSJed Brown           ierr = KSPSetDMActive(ilink->ksp,PETSC_FALSE);CHKERRQ(ierr);
2572fa5ba8aSJed Brown         }
2582fa5ba8aSJed Brown         ierr = PetscFree(dms);CHKERRQ(ierr);
259ea6813d4SMatthew G Knepley       } else if (dmcomplex) {
260ea6813d4SMatthew G Knepley         /* Define fields */
261*d2822287SMatthew G Knepley         PetscSection section, sectionGlobal;
262*d2822287SMatthew G Knepley         PetscInt    *fieldSizes, **fieldIndices;
263*d2822287SMatthew G Knepley         PetscInt     nF, f, pStart, pEnd, p;
264*d2822287SMatthew G Knepley 
265*d2822287SMatthew G Knepley         ierr = DMComplexGetDefaultSection(pc->dm, &section);CHKERRQ(ierr);
266*d2822287SMatthew G Knepley         ierr = DMComplexGetDefaultGlobalSection(pc->dm, &sectionGlobal);CHKERRQ(ierr);
267*d2822287SMatthew G Knepley         ierr = PetscSectionGetNumFields(section, &nF);CHKERRQ(ierr);
268*d2822287SMatthew G Knepley         ierr = PetscMalloc2(nF,PetscInt,&fieldSizes,nF,PetscInt *,&fieldIndices);CHKERRQ(ierr);
269*d2822287SMatthew G Knepley         ierr = PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);CHKERRQ(ierr);
270*d2822287SMatthew G Knepley         for(f = 0; f < nF; ++f) {
271*d2822287SMatthew G Knepley           fieldSizes[f] = 0;
272*d2822287SMatthew G Knepley         }
273*d2822287SMatthew G Knepley         for(p = pStart; p < pEnd; ++p) {
274*d2822287SMatthew G Knepley           PetscInt gdof;
275*d2822287SMatthew G Knepley 
276*d2822287SMatthew G Knepley           ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr);
277*d2822287SMatthew G Knepley           if (gdof > 0) {
278*d2822287SMatthew G Knepley             for(f = 0; f < nF; ++f) {
279*d2822287SMatthew G Knepley               PetscInt fcdof;
280*d2822287SMatthew G Knepley 
281*d2822287SMatthew G Knepley               ierr = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr);
282*d2822287SMatthew G Knepley               fieldSizes[f] += fcdof;
283*d2822287SMatthew G Knepley             }
284*d2822287SMatthew G Knepley           }
285*d2822287SMatthew G Knepley         }
286*d2822287SMatthew G Knepley         for(f = 0; f < nF; ++f) {
287*d2822287SMatthew G Knepley           ierr = PetscMalloc(fieldSizes[f] * sizeof(PetscInt), &fieldIndices[f]);CHKERRQ(ierr);
288*d2822287SMatthew G Knepley           fieldSizes[f] = 0;
289*d2822287SMatthew G Knepley         }
290*d2822287SMatthew G Knepley         for(p = pStart; p < pEnd; ++p) {
291*d2822287SMatthew G Knepley           PetscInt gdof, goff;
292*d2822287SMatthew G Knepley 
293*d2822287SMatthew G Knepley           ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr);
294*d2822287SMatthew G Knepley           if (gdof > 0) {
295*d2822287SMatthew G Knepley             ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
296*d2822287SMatthew G Knepley             for(f = 0; f < nF; ++f) {
297*d2822287SMatthew G Knepley               PetscInt fcdof, fc;
298*d2822287SMatthew G Knepley 
299*d2822287SMatthew G Knepley               ierr = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr);
300*d2822287SMatthew G Knepley               for(fc = 0; fc < fcdof; ++fc, ++fieldSizes[f]) {
301*d2822287SMatthew G Knepley                 fieldIndices[f][fieldSizes[f]] = goff++;
302*d2822287SMatthew G Knepley               }
303*d2822287SMatthew G Knepley             }
304*d2822287SMatthew G Knepley           }
305*d2822287SMatthew G Knepley         }
306*d2822287SMatthew G Knepley         for(f = 0; f < nF; ++f) {
307*d2822287SMatthew G Knepley           IS          field;
308*d2822287SMatthew G Knepley           const char *fieldname;
309*d2822287SMatthew G Knepley 
310*d2822287SMatthew G Knepley           ierr = PetscSectionGetFieldName(section, f, &fieldname);CHKERRQ(ierr);
311*d2822287SMatthew G Knepley           ierr = ISCreateGeneral(((PetscObject) pc)->comm, fieldSizes[f], fieldIndices[f], PETSC_OWN_POINTER, &field);CHKERRQ(ierr);
312*d2822287SMatthew G Knepley           ierr = PetscPrintf(((PetscObject) pc)->comm, "Field %s\n", fieldname);CHKERRQ(ierr);
313*d2822287SMatthew G Knepley           ierr = ISView(field, PETSC_NULL);CHKERRQ(ierr);
314*d2822287SMatthew G Knepley           ierr = PCFieldSplitSetIS(pc, fieldname, field);CHKERRQ(ierr);
315*d2822287SMatthew G Knepley           ierr = ISDestroy(&field);CHKERRQ(ierr);
316*d2822287SMatthew G Knepley         }
317*d2822287SMatthew G Knepley         ierr = PetscFree2(fieldSizes,fieldIndices);CHKERRQ(ierr);
3188b8307b2SJed Brown       }
31966ffff09SJed Brown     } else {
320521d7252SBarry Smith       if (jac->bs <= 0) {
321704ba839SBarry Smith         if (pc->pmat) {
322521d7252SBarry Smith           ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
323704ba839SBarry Smith         } else {
324704ba839SBarry Smith           jac->bs = 1;
325704ba839SBarry Smith         }
326521d7252SBarry Smith       }
327d32f9abdSBarry Smith 
328acfcf0e5SJed Brown       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
3296ce1633cSBarry Smith       if (stokes) {
3306ce1633cSBarry Smith         IS       zerodiags,rest;
3316ce1633cSBarry Smith         PetscInt nmin,nmax;
3326ce1633cSBarry Smith 
3336ce1633cSBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
3346ce1633cSBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
3356ce1633cSBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
3366ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
3376ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr);
338fcfd50ebSBarry Smith         ierr = ISDestroy(&zerodiags);CHKERRQ(ierr);
339fcfd50ebSBarry Smith         ierr = ISDestroy(&rest);CHKERRQ(ierr);
3406ce1633cSBarry Smith       } else {
341d32f9abdSBarry Smith         if (!flg) {
342d32f9abdSBarry Smith           /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
343d32f9abdSBarry Smith            then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
3446c924f48SJed Brown           ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
3456c924f48SJed Brown           if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
346d32f9abdSBarry Smith         }
3476c924f48SJed Brown         if (flg || !jac->splitdefined) {
348d32f9abdSBarry Smith           ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
349db4c96c1SJed Brown           for (i=0; i<jac->bs; i++) {
3506c924f48SJed Brown             char splitname[8];
3516c924f48SJed Brown             ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
352db4c96c1SJed Brown             ierr = PCFieldSplitSetFields(pc,splitname,1,&i);CHKERRQ(ierr);
35379416396SBarry Smith           }
35497bbdb24SBarry Smith           jac->defaultsplit = PETSC_TRUE;
355521d7252SBarry Smith         }
35666ffff09SJed Brown       }
3576ce1633cSBarry Smith     }
358edf189efSBarry Smith   } else if (jac->nsplits == 1) {
359edf189efSBarry Smith     if (ilink->is) {
360edf189efSBarry Smith       IS       is2;
361edf189efSBarry Smith       PetscInt nmin,nmax;
362edf189efSBarry Smith 
363edf189efSBarry Smith       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
364edf189efSBarry Smith       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
365db4c96c1SJed Brown       ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr);
366fcfd50ebSBarry Smith       ierr = ISDestroy(&is2);CHKERRQ(ierr);
367ea6813d4SMatthew G Knepley     } else SETERRQ(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
36863ec74ffSBarry Smith   } else if (jac->reset) {
36963ec74ffSBarry Smith     /* PCReset() has been called on this PC, ilink exists but all IS and Vec data structures in it must be rebuilt
370d0af7cd3SBarry Smith        This is basically the !ilink portion of code above copied from above and the allocation of the ilinks removed
371d0af7cd3SBarry Smith        since they already exist. This should be totally rewritten */
372d0af7cd3SBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
373d0af7cd3SBarry Smith     if (pc->dm && !stokes) {
374d0af7cd3SBarry Smith       PetscBool dmcomposite;
375d0af7cd3SBarry Smith       ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr);
376d0af7cd3SBarry Smith       if (dmcomposite) {
377d0af7cd3SBarry Smith         PetscInt nDM;
378d0af7cd3SBarry Smith         IS       *fields;
3792fa5ba8aSJed Brown         DM       *dms;
380d0af7cd3SBarry Smith         ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
381d0af7cd3SBarry Smith         ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr);
382d0af7cd3SBarry Smith         ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr);
3832fa5ba8aSJed Brown         ierr = PetscMalloc(nDM*sizeof(DM),&dms);CHKERRQ(ierr);
3842fa5ba8aSJed Brown         ierr = DMCompositeGetEntriesArray(pc->dm,dms);CHKERRQ(ierr);
385d0af7cd3SBarry Smith         for (i=0; i<nDM; i++) {
3862fa5ba8aSJed Brown           ierr = KSPSetDM(ilink->ksp,dms[i]);CHKERRQ(ierr);
3872fa5ba8aSJed Brown           ierr = KSPSetDMActive(ilink->ksp,PETSC_FALSE);CHKERRQ(ierr);
388d0af7cd3SBarry Smith           ilink->is = fields[i];
389d0af7cd3SBarry Smith           ilink     = ilink->next;
390edf189efSBarry Smith         }
391d0af7cd3SBarry Smith         ierr = PetscFree(fields);CHKERRQ(ierr);
3922fa5ba8aSJed Brown         ierr = PetscFree(dms);CHKERRQ(ierr);
393d0af7cd3SBarry Smith       }
394d0af7cd3SBarry Smith     } else {
395d0af7cd3SBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
396d0af7cd3SBarry Smith       if (stokes) {
397d0af7cd3SBarry Smith         IS       zerodiags,rest;
398d0af7cd3SBarry Smith         PetscInt nmin,nmax;
399d0af7cd3SBarry Smith 
400d0af7cd3SBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
401d0af7cd3SBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
402d0af7cd3SBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
40320b26d62SBarry Smith         ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
40420b26d62SBarry Smith         ierr = ISDestroy(&ilink->next->is);CHKERRQ(ierr);
405d0af7cd3SBarry Smith         ilink->is       = rest;
406d0af7cd3SBarry Smith         ilink->next->is = zerodiags;
40720b26d62SBarry Smith       } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used");
408d0af7cd3SBarry Smith     }
409d0af7cd3SBarry Smith   }
410d0af7cd3SBarry Smith 
411ea6813d4SMatthew 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);
41269a612a9SBarry Smith   PetscFunctionReturn(0);
41369a612a9SBarry Smith }
41469a612a9SBarry Smith 
41569a612a9SBarry Smith #undef __FUNCT__
41669a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
41769a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
41869a612a9SBarry Smith {
41969a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
42069a612a9SBarry Smith   PetscErrorCode    ierr;
4215a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
422cf502942SBarry Smith   PetscInt          i,nsplit,ccsize;
42369a612a9SBarry Smith   MatStructure      flag = pc->flag;
424ace3abfcSBarry Smith   PetscBool         sorted;
42569a612a9SBarry Smith 
42669a612a9SBarry Smith   PetscFunctionBegin;
42769a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
42897bbdb24SBarry Smith   nsplit = jac->nsplits;
4295a9f2f41SSatish Balay   ilink  = jac->head;
43097bbdb24SBarry Smith 
43197bbdb24SBarry Smith   /* get the matrices for each split */
432704ba839SBarry Smith   if (!jac->issetup) {
4331b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
43497bbdb24SBarry Smith 
435704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
436704ba839SBarry Smith 
437704ba839SBarry Smith     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
43851f519a2SBarry Smith     bs     = jac->bs;
43997bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
440cf502942SBarry Smith     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
4411b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
4421b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
4431b9fc7fcSBarry Smith       if (jac->defaultsplit) {
444704ba839SBarry Smith         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
445704ba839SBarry Smith       } else if (!ilink->is) {
446ccb205f8SBarry Smith         if (ilink->nfields > 1) {
4475a9f2f41SSatish Balay           PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
4485a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
4491b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
4501b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
4511b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
45297bbdb24SBarry Smith             }
45397bbdb24SBarry Smith           }
454d67e408aSBarry Smith           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr);
455ccb205f8SBarry Smith         } else {
456704ba839SBarry Smith           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
457ccb205f8SBarry Smith         }
4583e197d65SBarry Smith       }
459704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
460e32f2f54SBarry Smith       if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
461704ba839SBarry Smith       ilink = ilink->next;
4621b9fc7fcSBarry Smith     }
4631b9fc7fcSBarry Smith   }
4641b9fc7fcSBarry Smith 
465704ba839SBarry Smith   ilink  = jac->head;
46697bbdb24SBarry Smith   if (!jac->pmat) {
467cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
468cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
4694aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
470704ba839SBarry Smith       ilink = ilink->next;
471cf502942SBarry Smith     }
47297bbdb24SBarry Smith   } else {
473cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
4744aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
475704ba839SBarry Smith       ilink = ilink->next;
476cf502942SBarry Smith     }
47797bbdb24SBarry Smith   }
478519d70e2SJed Brown   if (jac->realdiagonal) {
479519d70e2SJed Brown     ilink = jac->head;
480519d70e2SJed Brown     if (!jac->mat) {
481519d70e2SJed Brown       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr);
482519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
483519d70e2SJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
484519d70e2SJed Brown         ilink = ilink->next;
485519d70e2SJed Brown       }
486519d70e2SJed Brown     } else {
487519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
488966e74d7SJed Brown         if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);}
489519d70e2SJed Brown         ilink = ilink->next;
490519d70e2SJed Brown       }
491519d70e2SJed Brown     }
492519d70e2SJed Brown   } else {
493519d70e2SJed Brown     jac->mat = jac->pmat;
494519d70e2SJed Brown   }
49597bbdb24SBarry Smith 
4966c8605c2SJed Brown   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
49768dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
49868dd23aaSBarry Smith     ilink  = jac->head;
49968dd23aaSBarry Smith     if (!jac->Afield) {
50068dd23aaSBarry Smith       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
50168dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
5024aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
50368dd23aaSBarry Smith         ilink = ilink->next;
50468dd23aaSBarry Smith       }
50568dd23aaSBarry Smith     } else {
50668dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
5074aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
50868dd23aaSBarry Smith         ilink = ilink->next;
50968dd23aaSBarry Smith       }
51068dd23aaSBarry Smith     }
51168dd23aaSBarry Smith   }
51268dd23aaSBarry Smith 
5133b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
5143b224e63SBarry Smith     IS       ccis;
5154aa3045dSJed Brown     PetscInt rstart,rend;
516e7e72b3dSBarry Smith     if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
51768dd23aaSBarry Smith 
518e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
519e6cab6aaSJed Brown     ierr  = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
520e6cab6aaSJed Brown 
5213b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
5223b224e63SBarry Smith     if (jac->schur) {
5233b224e63SBarry Smith       ilink = jac->head;
5244aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
5254aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
526fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
5273b224e63SBarry Smith       ilink = ilink->next;
5284aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
5294aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
530fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
531519d70e2SJed Brown       ierr  = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr);
532084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr);
5333b224e63SBarry Smith 
5343b224e63SBarry Smith      } else {
5351cee3971SBarry Smith       KSP ksp;
5366c924f48SJed Brown       char schurprefix[256];
5373b224e63SBarry Smith 
538a04f6461SBarry Smith       /* extract the A01 and A10 matrices */
5393b224e63SBarry Smith       ilink = jac->head;
5404aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
5414aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
542fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
5433b224e63SBarry Smith       ilink = ilink->next;
5444aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
5454aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
546fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
5477d50072eSJed Brown       /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */
548519d70e2SJed Brown       ierr  = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr);
549a04f6461SBarry Smith       /* set tabbing and options prefix of KSP inside the MatSchur */
5501cee3971SBarry Smith       ierr  = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
551e69d4d44SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr);
552a04f6461SBarry Smith       ierr  = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",jac->head->splitname);CHKERRQ(ierr);
553a04f6461SBarry Smith       ierr  = KSPSetOptionsPrefix(ksp,schurprefix);CHKERRQ(ierr);
55420b26d62SBarry Smith       /* Need to call this everytime because new matrix is being created */
55520b26d62SBarry Smith       ierr  = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
5567ae8954aSSatish Balay       ierr  = MatSetUp(jac->schur);CHKERRQ(ierr);
5573b224e63SBarry Smith 
5583b224e63SBarry Smith       ierr  = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr);
5599005cf84SBarry Smith       ierr  = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr);
5601cee3971SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
561084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
562084e4875SJed Brown       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
563084e4875SJed Brown         PC pc;
564cd5f4a64SJed Brown         ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr);
565084e4875SJed Brown         ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
566084e4875SJed Brown         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
567e69d4d44SBarry Smith       }
5686c924f48SJed Brown       ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
5696c924f48SJed Brown       ierr  = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr);
5703b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
57120b26d62SBarry Smith       /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */
57220b26d62SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
5733b224e63SBarry Smith 
5743b224e63SBarry Smith       ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr);
5753b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr);
5763b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr);
5773b224e63SBarry Smith       ilink = jac->head;
5783b224e63SBarry Smith       ilink->x = jac->x[0]; ilink->y = jac->y[0];
5793b224e63SBarry Smith       ilink = ilink->next;
5803b224e63SBarry Smith       ilink->x = jac->x[1]; ilink->y = jac->y[1];
5813b224e63SBarry Smith     }
5823b224e63SBarry Smith   } else {
58397bbdb24SBarry Smith     /* set up the individual PCs */
58497bbdb24SBarry Smith     i    = 0;
5855a9f2f41SSatish Balay     ilink = jac->head;
5865a9f2f41SSatish Balay     while (ilink) {
587519d70e2SJed Brown       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr);
5883b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
589c1570756SJed Brown       if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);}
5905a9f2f41SSatish Balay       ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
59197bbdb24SBarry Smith       i++;
5925a9f2f41SSatish Balay       ilink = ilink->next;
5930971522cSBarry Smith     }
59497bbdb24SBarry Smith 
59597bbdb24SBarry Smith     /* create work vectors for each split */
5961b9fc7fcSBarry Smith     if (!jac->x) {
59797bbdb24SBarry Smith       ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
5985a9f2f41SSatish Balay       ilink = jac->head;
59997bbdb24SBarry Smith       for (i=0; i<nsplit; i++) {
600906ed7ccSBarry Smith         Vec *vl,*vr;
601906ed7ccSBarry Smith 
602906ed7ccSBarry Smith         ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
603906ed7ccSBarry Smith         ilink->x  = *vr;
604906ed7ccSBarry Smith         ilink->y  = *vl;
605906ed7ccSBarry Smith         ierr      = PetscFree(vr);CHKERRQ(ierr);
606906ed7ccSBarry Smith         ierr      = PetscFree(vl);CHKERRQ(ierr);
6075a9f2f41SSatish Balay         jac->x[i] = ilink->x;
6085a9f2f41SSatish Balay         jac->y[i] = ilink->y;
6095a9f2f41SSatish Balay         ilink     = ilink->next;
61097bbdb24SBarry Smith       }
6113b224e63SBarry Smith     }
6123b224e63SBarry Smith   }
6133b224e63SBarry Smith 
6143b224e63SBarry Smith 
615d0f46423SBarry Smith   if (!jac->head->sctx) {
6163b224e63SBarry Smith     Vec xtmp;
6173b224e63SBarry Smith 
61879416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
6191b9fc7fcSBarry Smith 
6205a9f2f41SSatish Balay     ilink = jac->head;
6211b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
6221b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
623704ba839SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
6245a9f2f41SSatish Balay       ilink = ilink->next;
6251b9fc7fcSBarry Smith     }
6266bf464f9SBarry Smith     ierr = VecDestroy(&xtmp);CHKERRQ(ierr);
6271b9fc7fcSBarry Smith   }
628c1570756SJed Brown   jac->suboptionsset = PETSC_TRUE;
6290971522cSBarry Smith   PetscFunctionReturn(0);
6300971522cSBarry Smith }
6310971522cSBarry Smith 
6325a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
633ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
634ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
6355a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
636ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
637ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
63879416396SBarry Smith 
6390971522cSBarry Smith #undef __FUNCT__
6403b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur"
6413b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
6423b224e63SBarry Smith {
6433b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
6443b224e63SBarry Smith   PetscErrorCode    ierr;
6453b224e63SBarry Smith   KSP               ksp;
6463b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
6473b224e63SBarry Smith 
6483b224e63SBarry Smith   PetscFunctionBegin;
6493b224e63SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
6503b224e63SBarry Smith 
651c5d2311dSJed Brown   switch (jac->schurfactorization) {
652c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG:
653a04f6461SBarry Smith       /* [A00 0; 0 -S], positive definite, suitable for MINRES */
654c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
655c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
656c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
657c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
658c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
659c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
660c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
661c5d2311dSJed Brown       ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr);
662c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
663c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
664c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
665c5d2311dSJed Brown       break;
666c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER:
667a04f6461SBarry Smith       /* [A00 0; A10 S], suitable for left preconditioning */
668c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
669c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
670c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
671c5d2311dSJed Brown       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
672c5d2311dSJed Brown       ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr);
673c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
674c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
675c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
676c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
677c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
678c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
679c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
680c5d2311dSJed Brown       break;
681c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER:
682a04f6461SBarry Smith       /* [A00 A01; 0 S], suitable for right preconditioning */
683c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
684c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
685c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
686c5d2311dSJed Brown       ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
687c5d2311dSJed Brown       ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr);
688c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
689c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
690c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
691c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
692c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
693c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
694c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
695c5d2311dSJed Brown       break;
696c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL:
697a04f6461SBarry 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 */
6983b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6993b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7003b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
7013b224e63SBarry Smith       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
7023b224e63SBarry Smith       ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
7033b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7043b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7053b224e63SBarry Smith 
7063b224e63SBarry Smith       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
7073b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7083b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7093b224e63SBarry Smith 
7103b224e63SBarry Smith       ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
7113b224e63SBarry Smith       ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
7123b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
7133b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7143b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
715c5d2311dSJed Brown   }
7163b224e63SBarry Smith   PetscFunctionReturn(0);
7173b224e63SBarry Smith }
7183b224e63SBarry Smith 
7193b224e63SBarry Smith #undef __FUNCT__
7200971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
7210971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
7220971522cSBarry Smith {
7230971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
7240971522cSBarry Smith   PetscErrorCode    ierr;
7255a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
726af93645bSJed Brown   PetscInt          cnt;
7270971522cSBarry Smith 
7280971522cSBarry Smith   PetscFunctionBegin;
72951f519a2SBarry Smith   CHKMEMQ;
73051f519a2SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
73151f519a2SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
73251f519a2SBarry Smith 
73379416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
7341b9fc7fcSBarry Smith     if (jac->defaultsplit) {
7350971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
7365a9f2f41SSatish Balay       while (ilink) {
7375a9f2f41SSatish Balay         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
7385a9f2f41SSatish Balay         ilink = ilink->next;
7390971522cSBarry Smith       }
7400971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
7411b9fc7fcSBarry Smith     } else {
742efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
7435a9f2f41SSatish Balay       while (ilink) {
7445a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
7455a9f2f41SSatish Balay         ilink = ilink->next;
7461b9fc7fcSBarry Smith       }
7471b9fc7fcSBarry Smith     }
74816913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
74979416396SBarry Smith     if (!jac->w1) {
75079416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
75179416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
75279416396SBarry Smith     }
753efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
7545a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
7553e197d65SBarry Smith     cnt = 1;
7565a9f2f41SSatish Balay     while (ilink->next) {
7575a9f2f41SSatish Balay       ilink = ilink->next;
7583e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
7593e197d65SBarry Smith       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
7603e197d65SBarry Smith       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
7613e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7623e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7633e197d65SBarry Smith       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
7643e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7653e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
7663e197d65SBarry Smith     }
76751f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
76811755939SBarry Smith       cnt -= 2;
76951f519a2SBarry Smith       while (ilink->previous) {
77051f519a2SBarry Smith         ilink = ilink->previous;
77111755939SBarry Smith         /* compute the residual only over the part of the vector needed */
77211755939SBarry Smith         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
77311755939SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
77411755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
77511755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
77611755939SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
77711755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
77811755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
77951f519a2SBarry Smith       }
78011755939SBarry Smith     }
78165e19b50SBarry Smith   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
78251f519a2SBarry Smith   CHKMEMQ;
7830971522cSBarry Smith   PetscFunctionReturn(0);
7840971522cSBarry Smith }
7850971522cSBarry Smith 
786421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
787ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
788ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
789421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
790ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
791ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
792421e10b8SBarry Smith 
793421e10b8SBarry Smith #undef __FUNCT__
7948c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit"
795421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
796421e10b8SBarry Smith {
797421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
798421e10b8SBarry Smith   PetscErrorCode    ierr;
799421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
800421e10b8SBarry Smith 
801421e10b8SBarry Smith   PetscFunctionBegin;
802421e10b8SBarry Smith   CHKMEMQ;
803421e10b8SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
804421e10b8SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
805421e10b8SBarry Smith 
806421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
807421e10b8SBarry Smith     if (jac->defaultsplit) {
808421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
809421e10b8SBarry Smith       while (ilink) {
810421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
811421e10b8SBarry Smith 	ilink = ilink->next;
812421e10b8SBarry Smith       }
813421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
814421e10b8SBarry Smith     } else {
815421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
816421e10b8SBarry Smith       while (ilink) {
817421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
818421e10b8SBarry Smith 	ilink = ilink->next;
819421e10b8SBarry Smith       }
820421e10b8SBarry Smith     }
821421e10b8SBarry Smith   } else {
822421e10b8SBarry Smith     if (!jac->w1) {
823421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
824421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
825421e10b8SBarry Smith     }
826421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
827421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
828421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
829421e10b8SBarry Smith       while (ilink->next) {
830421e10b8SBarry Smith         ilink = ilink->next;
8319989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
832421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
833421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
834421e10b8SBarry Smith       }
835421e10b8SBarry Smith       while (ilink->previous) {
836421e10b8SBarry Smith         ilink = ilink->previous;
8379989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
838421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
839421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
840421e10b8SBarry Smith       }
841421e10b8SBarry Smith     } else {
842421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
843421e10b8SBarry Smith 	ilink = ilink->next;
844421e10b8SBarry Smith       }
845421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
846421e10b8SBarry Smith       while (ilink->previous) {
847421e10b8SBarry Smith 	ilink = ilink->previous;
8489989ab13SBarry Smith 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
849421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
850421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
851421e10b8SBarry Smith       }
852421e10b8SBarry Smith     }
853421e10b8SBarry Smith   }
854421e10b8SBarry Smith   CHKMEMQ;
855421e10b8SBarry Smith   PetscFunctionReturn(0);
856421e10b8SBarry Smith }
857421e10b8SBarry Smith 
8580971522cSBarry Smith #undef __FUNCT__
859574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit"
860574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc)
8610971522cSBarry Smith {
8620971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
8630971522cSBarry Smith   PetscErrorCode    ierr;
8645a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
8650971522cSBarry Smith 
8660971522cSBarry Smith   PetscFunctionBegin;
8675a9f2f41SSatish Balay   while (ilink) {
868574deadeSBarry Smith     ierr = KSPReset(ilink->ksp);CHKERRQ(ierr);
869fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->x);CHKERRQ(ierr);
870fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->y);CHKERRQ(ierr);
871fcfd50ebSBarry Smith     ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr);
872fcfd50ebSBarry Smith     ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
8735a9f2f41SSatish Balay     next = ilink->next;
8745a9f2f41SSatish Balay     ilink = next;
8750971522cSBarry Smith   }
87605b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
877574deadeSBarry Smith   if (jac->mat && jac->mat != jac->pmat) {
878574deadeSBarry Smith     ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);
879574deadeSBarry Smith   } else if (jac->mat) {
880574deadeSBarry Smith     jac->mat = PETSC_NULL;
881574deadeSBarry Smith   }
88297bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
88368dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
8846bf464f9SBarry Smith   ierr = VecDestroy(&jac->w1);CHKERRQ(ierr);
8856bf464f9SBarry Smith   ierr = VecDestroy(&jac->w2);CHKERRQ(ierr);
8866bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur);CHKERRQ(ierr);
8876bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
8886bf464f9SBarry Smith   ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr);
8896bf464f9SBarry Smith   ierr = MatDestroy(&jac->B);CHKERRQ(ierr);
8906bf464f9SBarry Smith   ierr = MatDestroy(&jac->C);CHKERRQ(ierr);
89163ec74ffSBarry Smith   jac->reset = PETSC_TRUE;
892574deadeSBarry Smith   PetscFunctionReturn(0);
893574deadeSBarry Smith }
894574deadeSBarry Smith 
895574deadeSBarry Smith #undef __FUNCT__
896574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
897574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
898574deadeSBarry Smith {
899574deadeSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
900574deadeSBarry Smith   PetscErrorCode    ierr;
901574deadeSBarry Smith   PC_FieldSplitLink ilink = jac->head,next;
902574deadeSBarry Smith 
903574deadeSBarry Smith   PetscFunctionBegin;
904574deadeSBarry Smith   ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr);
905574deadeSBarry Smith   while (ilink) {
9066bf464f9SBarry Smith     ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr);
907574deadeSBarry Smith     next = ilink->next;
908574deadeSBarry Smith     ierr = PetscFree(ilink->splitname);CHKERRQ(ierr);
909574deadeSBarry Smith     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
910574deadeSBarry Smith     ierr = PetscFree(ilink);CHKERRQ(ierr);
911574deadeSBarry Smith     ilink = next;
912574deadeSBarry Smith   }
913574deadeSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
914c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
915d88c663fSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","",PETSC_NULL);CHKERRQ(ierr);
916d88c663fSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","",PETSC_NULL);CHKERRQ(ierr);
917d88c663fSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","",PETSC_NULL);CHKERRQ(ierr);
918d88c663fSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","",PETSC_NULL);CHKERRQ(ierr);
919d88c663fSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","",PETSC_NULL);CHKERRQ(ierr);
9200971522cSBarry Smith   PetscFunctionReturn(0);
9210971522cSBarry Smith }
9220971522cSBarry Smith 
9230971522cSBarry Smith #undef __FUNCT__
9240971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
9250971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
9260971522cSBarry Smith {
9271b9fc7fcSBarry Smith   PetscErrorCode  ierr;
9286c924f48SJed Brown   PetscInt        bs;
929bc59fbc5SBarry Smith   PetscBool       flg,stokes = PETSC_FALSE;
9309dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
9313b224e63SBarry Smith   PCCompositeType ctype;
9321b9fc7fcSBarry Smith 
9330971522cSBarry Smith   PetscFunctionBegin;
9341b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
935acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr);
93651f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
93751f519a2SBarry Smith   if (flg) {
93851f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
93951f519a2SBarry Smith   }
940704ba839SBarry Smith 
941435f959eSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
942c0adfefeSBarry Smith   if (stokes) {
943c0adfefeSBarry Smith     ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
944c0adfefeSBarry Smith     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
945c0adfefeSBarry Smith   }
946c0adfefeSBarry Smith 
9473b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
9483b224e63SBarry Smith   if (flg) {
9493b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
9503b224e63SBarry Smith   }
951d32f9abdSBarry Smith 
952c30613efSMatthew Knepley   /* Only setup fields once */
953c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
954d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
955d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
9566c924f48SJed Brown     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
9576c924f48SJed Brown     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
958d32f9abdSBarry Smith   }
959c5d2311dSJed Brown   if (jac->type == PC_COMPOSITE_SCHUR) {
960c5d2311dSJed Brown     ierr = PetscOptionsEnum("-pc_fieldsplit_schur_factorization_type","Factorization to use","None",PCFieldSplitSchurFactorizationTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,PETSC_NULL);CHKERRQ(ierr);
961084e4875SJed 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);
962c5d2311dSJed Brown   }
9631b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
9640971522cSBarry Smith   PetscFunctionReturn(0);
9650971522cSBarry Smith }
9660971522cSBarry Smith 
9670971522cSBarry Smith /*------------------------------------------------------------------------------------*/
9680971522cSBarry Smith 
9690971522cSBarry Smith EXTERN_C_BEGIN
9700971522cSBarry Smith #undef __FUNCT__
9710971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
9727087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields)
9730971522cSBarry Smith {
97497bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
9750971522cSBarry Smith   PetscErrorCode    ierr;
9765a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
97769a612a9SBarry Smith   char              prefix[128];
97851f519a2SBarry Smith   PetscInt          i;
9790971522cSBarry Smith 
9800971522cSBarry Smith   PetscFunctionBegin;
9816c924f48SJed Brown   if (jac->splitdefined) {
9826c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
9836c924f48SJed Brown     PetscFunctionReturn(0);
9846c924f48SJed Brown   }
98551f519a2SBarry Smith   for (i=0; i<n; i++) {
986e32f2f54SBarry 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);
987e32f2f54SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
98851f519a2SBarry Smith   }
989704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
990a04f6461SBarry Smith   if (splitname) {
991db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
992a04f6461SBarry Smith   } else {
993a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
994a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
995a04f6461SBarry Smith   }
996704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
9975a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
9985a9f2f41SSatish Balay   ilink->nfields = n;
9995a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
10007adad957SLisandro Dalcin   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
10011cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
10025a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
10039005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
100469a612a9SBarry Smith 
1005a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
10065a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
10070971522cSBarry Smith 
10080971522cSBarry Smith   if (!next) {
10095a9f2f41SSatish Balay     jac->head       = ilink;
101051f519a2SBarry Smith     ilink->previous = PETSC_NULL;
10110971522cSBarry Smith   } else {
10120971522cSBarry Smith     while (next->next) {
10130971522cSBarry Smith       next = next->next;
10140971522cSBarry Smith     }
10155a9f2f41SSatish Balay     next->next      = ilink;
101651f519a2SBarry Smith     ilink->previous = next;
10170971522cSBarry Smith   }
10180971522cSBarry Smith   jac->nsplits++;
10190971522cSBarry Smith   PetscFunctionReturn(0);
10200971522cSBarry Smith }
10210971522cSBarry Smith EXTERN_C_END
10220971522cSBarry Smith 
1023e69d4d44SBarry Smith EXTERN_C_BEGIN
1024e69d4d44SBarry Smith #undef __FUNCT__
1025e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
10267087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
1027e69d4d44SBarry Smith {
1028e69d4d44SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
1029e69d4d44SBarry Smith   PetscErrorCode ierr;
1030e69d4d44SBarry Smith 
1031e69d4d44SBarry Smith   PetscFunctionBegin;
1032519d70e2SJed Brown   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
1033e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
1034e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
103513e0d083SBarry Smith   if (n) *n = jac->nsplits;
1036e69d4d44SBarry Smith   PetscFunctionReturn(0);
1037e69d4d44SBarry Smith }
1038e69d4d44SBarry Smith EXTERN_C_END
10390971522cSBarry Smith 
10400971522cSBarry Smith EXTERN_C_BEGIN
10410971522cSBarry Smith #undef __FUNCT__
104269a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
10437087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
10440971522cSBarry Smith {
10450971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
10460971522cSBarry Smith   PetscErrorCode    ierr;
10470971522cSBarry Smith   PetscInt          cnt = 0;
10485a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
10490971522cSBarry Smith 
10500971522cSBarry Smith   PetscFunctionBegin;
10515d480477SMatthew G Knepley   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
10525a9f2f41SSatish Balay   while (ilink) {
10535a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
10545a9f2f41SSatish Balay     ilink = ilink->next;
10550971522cSBarry Smith   }
10565d480477SMatthew 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);
105713e0d083SBarry Smith   if (n) *n = jac->nsplits;
10580971522cSBarry Smith   PetscFunctionReturn(0);
10590971522cSBarry Smith }
10600971522cSBarry Smith EXTERN_C_END
10610971522cSBarry Smith 
1062704ba839SBarry Smith EXTERN_C_BEGIN
1063704ba839SBarry Smith #undef __FUNCT__
1064704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
10657087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
1066704ba839SBarry Smith {
1067704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1068704ba839SBarry Smith   PetscErrorCode    ierr;
1069704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
1070704ba839SBarry Smith   char              prefix[128];
1071704ba839SBarry Smith 
1072704ba839SBarry Smith   PetscFunctionBegin;
10736c924f48SJed Brown   if (jac->splitdefined) {
10746c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
10756c924f48SJed Brown     PetscFunctionReturn(0);
10766c924f48SJed Brown   }
107716913363SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
1078a04f6461SBarry Smith   if (splitname) {
1079db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1080a04f6461SBarry Smith   } else {
1081a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
1082d88c663fSJed Brown     ierr = PetscSNPrintf(ilink->splitname,3,"%D",jac->nsplits);CHKERRQ(ierr);
1083a04f6461SBarry Smith   }
10841ab39975SBarry Smith   ilink->is      = is;
10851ab39975SBarry Smith   ierr           = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1086704ba839SBarry Smith   ilink->next    = PETSC_NULL;
1087704ba839SBarry Smith   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
10881cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1089704ba839SBarry Smith   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
10909005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
1091704ba839SBarry Smith 
1092a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
1093704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
1094704ba839SBarry Smith 
1095704ba839SBarry Smith   if (!next) {
1096704ba839SBarry Smith     jac->head       = ilink;
1097704ba839SBarry Smith     ilink->previous = PETSC_NULL;
1098704ba839SBarry Smith   } else {
1099704ba839SBarry Smith     while (next->next) {
1100704ba839SBarry Smith       next = next->next;
1101704ba839SBarry Smith     }
1102704ba839SBarry Smith     next->next      = ilink;
1103704ba839SBarry Smith     ilink->previous = next;
1104704ba839SBarry Smith   }
1105704ba839SBarry Smith   jac->nsplits++;
1106704ba839SBarry Smith 
1107704ba839SBarry Smith   PetscFunctionReturn(0);
1108704ba839SBarry Smith }
1109704ba839SBarry Smith EXTERN_C_END
1110704ba839SBarry Smith 
11110971522cSBarry Smith #undef __FUNCT__
11120971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
11130971522cSBarry Smith /*@
11140971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
11150971522cSBarry Smith 
1116ad4df100SBarry Smith     Logically Collective on PC
11170971522cSBarry Smith 
11180971522cSBarry Smith     Input Parameters:
11190971522cSBarry Smith +   pc  - the preconditioner context
1120a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
11210971522cSBarry Smith .   n - the number of fields in this split
1122db4c96c1SJed Brown -   fields - the fields in this split
11230971522cSBarry Smith 
11240971522cSBarry Smith     Level: intermediate
11250971522cSBarry Smith 
1126d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
1127d32f9abdSBarry Smith 
1128d32f9abdSBarry Smith      The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block
1129d32f9abdSBarry 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
1130d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1131d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
1132d32f9abdSBarry Smith 
1133db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
1134db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1135db4c96c1SJed Brown 
1136d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
11370971522cSBarry Smith 
11380971522cSBarry Smith @*/
11397087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields)
11400971522cSBarry Smith {
11414ac538c5SBarry Smith   PetscErrorCode ierr;
11420971522cSBarry Smith 
11430971522cSBarry Smith   PetscFunctionBegin;
11440700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1145db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1146db4c96c1SJed Brown   if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1147db4c96c1SJed Brown   PetscValidIntPointer(fields,3);
11484ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *),(pc,splitname,n,fields));CHKERRQ(ierr);
11490971522cSBarry Smith   PetscFunctionReturn(0);
11500971522cSBarry Smith }
11510971522cSBarry Smith 
11520971522cSBarry Smith #undef __FUNCT__
1153704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
1154704ba839SBarry Smith /*@
1155704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
1156704ba839SBarry Smith 
1157ad4df100SBarry Smith     Logically Collective on PC
1158704ba839SBarry Smith 
1159704ba839SBarry Smith     Input Parameters:
1160704ba839SBarry Smith +   pc  - the preconditioner context
1161a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
1162db4c96c1SJed Brown -   is - the index set that defines the vector elements in this field
1163704ba839SBarry Smith 
1164d32f9abdSBarry Smith 
1165a6ffb8dbSJed Brown     Notes:
1166a6ffb8dbSJed Brown     Use PCFieldSplitSetFields(), for fields defined by strided types.
1167a6ffb8dbSJed Brown 
1168db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
1169db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1170d32f9abdSBarry Smith 
1171704ba839SBarry Smith     Level: intermediate
1172704ba839SBarry Smith 
1173704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
1174704ba839SBarry Smith 
1175704ba839SBarry Smith @*/
11767087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1177704ba839SBarry Smith {
11784ac538c5SBarry Smith   PetscErrorCode ierr;
1179704ba839SBarry Smith 
1180704ba839SBarry Smith   PetscFunctionBegin;
11810700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1182d88c663fSJed Brown   if (splitname) PetscValidCharPointer(splitname,2);
1183db4c96c1SJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,3);
11844ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr);
1185704ba839SBarry Smith   PetscFunctionReturn(0);
1186704ba839SBarry Smith }
1187704ba839SBarry Smith 
1188704ba839SBarry Smith #undef __FUNCT__
118957a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS"
119057a9adfeSMatthew G Knepley /*@
119157a9adfeSMatthew G Knepley     PCFieldSplitGetIS - Retrieves the elements for a field as an IS
119257a9adfeSMatthew G Knepley 
119357a9adfeSMatthew G Knepley     Logically Collective on PC
119457a9adfeSMatthew G Knepley 
119557a9adfeSMatthew G Knepley     Input Parameters:
119657a9adfeSMatthew G Knepley +   pc  - the preconditioner context
119757a9adfeSMatthew G Knepley -   splitname - name of this split
119857a9adfeSMatthew G Knepley 
119957a9adfeSMatthew G Knepley     Output Parameter:
120057a9adfeSMatthew G Knepley -   is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found
120157a9adfeSMatthew G Knepley 
120257a9adfeSMatthew G Knepley     Level: intermediate
120357a9adfeSMatthew G Knepley 
120457a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS()
120557a9adfeSMatthew G Knepley 
120657a9adfeSMatthew G Knepley @*/
120757a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
120857a9adfeSMatthew G Knepley {
120957a9adfeSMatthew G Knepley   PetscErrorCode ierr;
121057a9adfeSMatthew G Knepley 
121157a9adfeSMatthew G Knepley   PetscFunctionBegin;
121257a9adfeSMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
121357a9adfeSMatthew G Knepley   PetscValidCharPointer(splitname,2);
121457a9adfeSMatthew G Knepley   PetscValidPointer(is,3);
121557a9adfeSMatthew G Knepley   {
121657a9adfeSMatthew G Knepley     PC_FieldSplit    *jac   = (PC_FieldSplit *) pc->data;
121757a9adfeSMatthew G Knepley     PC_FieldSplitLink ilink = jac->head;
121857a9adfeSMatthew G Knepley     PetscBool         found;
121957a9adfeSMatthew G Knepley 
122057a9adfeSMatthew G Knepley     *is = PETSC_NULL;
122157a9adfeSMatthew G Knepley     while(ilink) {
122257a9adfeSMatthew G Knepley       ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr);
122357a9adfeSMatthew G Knepley       if (found) {
122457a9adfeSMatthew G Knepley         *is = ilink->is;
122557a9adfeSMatthew G Knepley         break;
122657a9adfeSMatthew G Knepley       }
122757a9adfeSMatthew G Knepley       ilink = ilink->next;
122857a9adfeSMatthew G Knepley     }
122957a9adfeSMatthew G Knepley   }
123057a9adfeSMatthew G Knepley   PetscFunctionReturn(0);
123157a9adfeSMatthew G Knepley }
123257a9adfeSMatthew G Knepley 
123357a9adfeSMatthew G Knepley #undef __FUNCT__
123451f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
123551f519a2SBarry Smith /*@
123651f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
123751f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
123851f519a2SBarry Smith 
1239ad4df100SBarry Smith     Logically Collective on PC
124051f519a2SBarry Smith 
124151f519a2SBarry Smith     Input Parameters:
124251f519a2SBarry Smith +   pc  - the preconditioner context
124351f519a2SBarry Smith -   bs - the block size
124451f519a2SBarry Smith 
124551f519a2SBarry Smith     Level: intermediate
124651f519a2SBarry Smith 
124751f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
124851f519a2SBarry Smith 
124951f519a2SBarry Smith @*/
12507087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
125151f519a2SBarry Smith {
12524ac538c5SBarry Smith   PetscErrorCode ierr;
125351f519a2SBarry Smith 
125451f519a2SBarry Smith   PetscFunctionBegin;
12550700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1256c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,bs,2);
12574ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr);
125851f519a2SBarry Smith   PetscFunctionReturn(0);
125951f519a2SBarry Smith }
126051f519a2SBarry Smith 
126151f519a2SBarry Smith #undef __FUNCT__
126269a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
12630971522cSBarry Smith /*@C
126469a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
12650971522cSBarry Smith 
126669a612a9SBarry Smith    Collective on KSP
12670971522cSBarry Smith 
12680971522cSBarry Smith    Input Parameter:
12690971522cSBarry Smith .  pc - the preconditioner context
12700971522cSBarry Smith 
12710971522cSBarry Smith    Output Parameters:
127213e0d083SBarry Smith +  n - the number of splits
127369a612a9SBarry Smith -  pc - the array of KSP contexts
12740971522cSBarry Smith 
12750971522cSBarry Smith    Note:
1276d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
1277d32f9abdSBarry Smith    (not the KSP just the array that contains them).
12780971522cSBarry Smith 
127969a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
12800971522cSBarry Smith 
12810971522cSBarry Smith    Level: advanced
12820971522cSBarry Smith 
12830971522cSBarry Smith .seealso: PCFIELDSPLIT
12840971522cSBarry Smith @*/
12857087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
12860971522cSBarry Smith {
12874ac538c5SBarry Smith   PetscErrorCode ierr;
12880971522cSBarry Smith 
12890971522cSBarry Smith   PetscFunctionBegin;
12900700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
129113e0d083SBarry Smith   if (n) PetscValidIntPointer(n,2);
12924ac538c5SBarry Smith   ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
12930971522cSBarry Smith   PetscFunctionReturn(0);
12940971522cSBarry Smith }
12950971522cSBarry Smith 
1296e69d4d44SBarry Smith #undef __FUNCT__
1297e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition"
1298e69d4d44SBarry Smith /*@
1299e69d4d44SBarry Smith     PCFieldSplitSchurPrecondition -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
1300a04f6461SBarry Smith       A11 matrix. Otherwise no preconditioner is used.
1301e69d4d44SBarry Smith 
1302e69d4d44SBarry Smith     Collective on PC
1303e69d4d44SBarry Smith 
1304e69d4d44SBarry Smith     Input Parameters:
1305e69d4d44SBarry Smith +   pc  - the preconditioner context
1306fd1303e9SJungho Lee .   ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default
1307084e4875SJed Brown -   userpre - matrix to use for preconditioning, or PETSC_NULL
1308084e4875SJed Brown 
1309e69d4d44SBarry Smith     Options Database:
1310084e4875SJed Brown .     -pc_fieldsplit_schur_precondition <self,user,diag> default is diag
1311e69d4d44SBarry Smith 
1312fd1303e9SJungho Lee     Notes:
1313fd1303e9SJungho Lee $    If ptype is
1314fd1303e9SJungho Lee $        user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument
1315fd1303e9SJungho Lee $             to this function).
1316fd1303e9SJungho Lee $        diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original
1317fd1303e9SJungho Lee $             matrix associated with the Schur complement (i.e. A11)
1318fd1303e9SJungho Lee $        self the preconditioner for the Schur complement is generated from the Schur complement matrix itself:
1319fd1303e9SJungho Lee $             The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC
1320fd1303e9SJungho Lee $             preconditioner
1321fd1303e9SJungho Lee 
1322fd1303e9SJungho Lee      When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense
1323fd1303e9SJungho Lee     with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
1324fd1303e9SJungho Lee     -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement.
1325fd1303e9SJungho Lee 
1326fd1303e9SJungho Lee     Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in
1327fd1303e9SJungho Lee      the name since it sets a proceedure to use.
1328fd1303e9SJungho Lee 
1329e69d4d44SBarry Smith     Level: intermediate
1330e69d4d44SBarry Smith 
1331fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC
1332e69d4d44SBarry Smith 
1333e69d4d44SBarry Smith @*/
13347087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1335e69d4d44SBarry Smith {
13364ac538c5SBarry Smith   PetscErrorCode ierr;
1337e69d4d44SBarry Smith 
1338e69d4d44SBarry Smith   PetscFunctionBegin;
13390700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13404ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr);
1341e69d4d44SBarry Smith   PetscFunctionReturn(0);
1342e69d4d44SBarry Smith }
1343e69d4d44SBarry Smith 
1344e69d4d44SBarry Smith EXTERN_C_BEGIN
1345e69d4d44SBarry Smith #undef __FUNCT__
1346e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
13477087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1348e69d4d44SBarry Smith {
1349e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1350084e4875SJed Brown   PetscErrorCode  ierr;
1351e69d4d44SBarry Smith 
1352e69d4d44SBarry Smith   PetscFunctionBegin;
1353084e4875SJed Brown   jac->schurpre = ptype;
1354084e4875SJed Brown   if (pre) {
13556bf464f9SBarry Smith     ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
1356084e4875SJed Brown     jac->schur_user = pre;
1357084e4875SJed Brown     ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1358084e4875SJed Brown   }
1359e69d4d44SBarry Smith   PetscFunctionReturn(0);
1360e69d4d44SBarry Smith }
1361e69d4d44SBarry Smith EXTERN_C_END
1362e69d4d44SBarry Smith 
136330ad9308SMatthew Knepley #undef __FUNCT__
136430ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
136530ad9308SMatthew Knepley /*@C
13668c03b21aSDmitry Karpeev    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
136730ad9308SMatthew Knepley 
136830ad9308SMatthew Knepley    Collective on KSP
136930ad9308SMatthew Knepley 
137030ad9308SMatthew Knepley    Input Parameter:
137130ad9308SMatthew Knepley .  pc - the preconditioner context
137230ad9308SMatthew Knepley 
137330ad9308SMatthew Knepley    Output Parameters:
1374a04f6461SBarry Smith +  A00 - the (0,0) block
1375a04f6461SBarry Smith .  A01 - the (0,1) block
1376a04f6461SBarry Smith .  A10 - the (1,0) block
1377a04f6461SBarry Smith -  A11 - the (1,1) block
137830ad9308SMatthew Knepley 
137930ad9308SMatthew Knepley    Level: advanced
138030ad9308SMatthew Knepley 
138130ad9308SMatthew Knepley .seealso: PCFIELDSPLIT
138230ad9308SMatthew Knepley @*/
1383a04f6461SBarry Smith PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
138430ad9308SMatthew Knepley {
138530ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit *) pc->data;
138630ad9308SMatthew Knepley 
138730ad9308SMatthew Knepley   PetscFunctionBegin;
13880700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1389c1235816SBarry Smith   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
1390a04f6461SBarry Smith   if (A00) *A00 = jac->pmat[0];
1391a04f6461SBarry Smith   if (A01) *A01 = jac->B;
1392a04f6461SBarry Smith   if (A10) *A10 = jac->C;
1393a04f6461SBarry Smith   if (A11) *A11 = jac->pmat[1];
139430ad9308SMatthew Knepley   PetscFunctionReturn(0);
139530ad9308SMatthew Knepley }
139630ad9308SMatthew Knepley 
139779416396SBarry Smith EXTERN_C_BEGIN
139879416396SBarry Smith #undef __FUNCT__
139979416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
14007087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
140179416396SBarry Smith {
140279416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1403e69d4d44SBarry Smith   PetscErrorCode ierr;
140479416396SBarry Smith 
140579416396SBarry Smith   PetscFunctionBegin;
140679416396SBarry Smith   jac->type = type;
14073b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
14083b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
14093b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
1410e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1411e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1412e69d4d44SBarry Smith 
14133b224e63SBarry Smith   } else {
14143b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
14153b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
1416e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
14179e7d6b0aSBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr);
14183b224e63SBarry Smith   }
141979416396SBarry Smith   PetscFunctionReturn(0);
142079416396SBarry Smith }
142179416396SBarry Smith EXTERN_C_END
142279416396SBarry Smith 
142351f519a2SBarry Smith EXTERN_C_BEGIN
142451f519a2SBarry Smith #undef __FUNCT__
142551f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
14267087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
142751f519a2SBarry Smith {
142851f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
142951f519a2SBarry Smith 
143051f519a2SBarry Smith   PetscFunctionBegin;
143165e19b50SBarry Smith   if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
143265e19b50SBarry 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);
143351f519a2SBarry Smith   jac->bs = bs;
143451f519a2SBarry Smith   PetscFunctionReturn(0);
143551f519a2SBarry Smith }
143651f519a2SBarry Smith EXTERN_C_END
143751f519a2SBarry Smith 
143879416396SBarry Smith #undef __FUNCT__
143979416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
1440bc08b0f1SBarry Smith /*@
144179416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
144279416396SBarry Smith 
144379416396SBarry Smith    Collective on PC
144479416396SBarry Smith 
144579416396SBarry Smith    Input Parameter:
144679416396SBarry Smith .  pc - the preconditioner context
144781540f2fSBarry Smith .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
144879416396SBarry Smith 
144979416396SBarry Smith    Options Database Key:
1450a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
145179416396SBarry Smith 
145279416396SBarry Smith    Level: Developer
145379416396SBarry Smith 
145479416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
145579416396SBarry Smith 
145679416396SBarry Smith .seealso: PCCompositeSetType()
145779416396SBarry Smith 
145879416396SBarry Smith @*/
14597087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
146079416396SBarry Smith {
14614ac538c5SBarry Smith   PetscErrorCode ierr;
146279416396SBarry Smith 
146379416396SBarry Smith   PetscFunctionBegin;
14640700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
14654ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr);
146679416396SBarry Smith   PetscFunctionReturn(0);
146779416396SBarry Smith }
146879416396SBarry Smith 
14690971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
14700971522cSBarry Smith /*MC
1471a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
1472a04f6461SBarry Smith                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.
14730971522cSBarry Smith 
1474edf189efSBarry Smith      To set options on the solvers for each block append -fieldsplit_ to all the PC
1475edf189efSBarry Smith         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
14760971522cSBarry Smith 
1477a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
147869a612a9SBarry Smith          and set the options directly on the resulting KSP object
14790971522cSBarry Smith 
14800971522cSBarry Smith    Level: intermediate
14810971522cSBarry Smith 
148279416396SBarry Smith    Options Database Keys:
148381540f2fSBarry Smith +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
148481540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
148581540f2fSBarry Smith                               been supplied explicitly by -pc_fieldsplit_%d_fields
148681540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
148781540f2fSBarry Smith .   -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative>
1488e69d4d44SBarry Smith .   -pc_fieldsplit_schur_precondition <true,false> default is true
1489435f959eSBarry Smith .   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver
149043890ad2SJed 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)
149143890ad2SJed Brown -   -fieldsplit_NAME_pc_* - control inner preconditioner, NAME has same semantics as above
149279416396SBarry Smith 
149343890ad2SJed Brown    Notes:
149443890ad2SJed Brown     Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1495d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
1496d32f9abdSBarry Smith 
1497d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
1498d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1499d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
1500d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
1501d32f9abdSBarry Smith 
1502a04f6461SBarry Smith      For the Schur complement preconditioner if J = ( A00 A01 )
1503a04f6461SBarry Smith                                                     ( A10 A11 )
1504e69d4d44SBarry Smith      the preconditioner is
1505a04f6461SBarry Smith               (I   -B ksp(A00)) ( inv(A00)   0  (I         0  )
1506a04f6461SBarry Smith               (0    I       ) (   0    ksp(S) ) (-A10 ksp(A00) I  )
1507a04f6461SBarry Smith      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of
1508a04f6461SBarry 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).
1509a04f6461SBarry Smith      For PCFieldSplitGetKSP() when field number is
1510edf189efSBarry Smith      0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
1511a04f6461SBarry Smith      A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
15127e8cb189SBarry Smith      option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc
1513e69d4d44SBarry Smith 
1514edf189efSBarry Smith      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
1515edf189efSBarry Smith      is used automatically for a second block.
1516edf189efSBarry Smith 
1517ff218e97SBarry Smith      The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1.
1518ff218e97SBarry Smith      Generally it should be used with the AIJ format.
1519ff218e97SBarry Smith 
1520ff218e97SBarry Smith      The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see,
1521ff218e97SBarry Smith      for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT
1522ff218e97SBarry Smith      inside a smoother resulting in "Distributive Smoothers".
15230716a85fSBarry Smith 
1524a541d17aSBarry Smith    Concepts: physics based preconditioners, block preconditioners
15250971522cSBarry Smith 
15267e8cb189SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
1527e69d4d44SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
15280971522cSBarry Smith M*/
15290971522cSBarry Smith 
15300971522cSBarry Smith EXTERN_C_BEGIN
15310971522cSBarry Smith #undef __FUNCT__
15320971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
15337087cfbeSBarry Smith PetscErrorCode  PCCreate_FieldSplit(PC pc)
15340971522cSBarry Smith {
15350971522cSBarry Smith   PetscErrorCode ierr;
15360971522cSBarry Smith   PC_FieldSplit  *jac;
15370971522cSBarry Smith 
15380971522cSBarry Smith   PetscFunctionBegin;
153938f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
15400971522cSBarry Smith   jac->bs        = -1;
15410971522cSBarry Smith   jac->nsplits   = 0;
15423e197d65SBarry Smith   jac->type      = PC_COMPOSITE_MULTIPLICATIVE;
1543e6cab6aaSJed Brown   jac->schurpre  = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
1544c5d2311dSJed Brown   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL;
154551f519a2SBarry Smith 
15460971522cSBarry Smith   pc->data     = (void*)jac;
15470971522cSBarry Smith 
15480971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
1549421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
15500971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
1551574deadeSBarry Smith   pc->ops->reset             = PCReset_FieldSplit;
15520971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
15530971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
15540971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
15550971522cSBarry Smith   pc->ops->applyrichardson   = 0;
15560971522cSBarry Smith 
155769a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
155869a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
15590971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
15600971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1561704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
1562704ba839SBarry Smith                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
156379416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
156479416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
156551f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
156651f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
15670971522cSBarry Smith   PetscFunctionReturn(0);
15680971522cSBarry Smith }
15690971522cSBarry Smith EXTERN_C_END
15700971522cSBarry Smith 
15710971522cSBarry Smith 
1572a541d17aSBarry Smith 
1573