xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision ff218e97a57ed641f3ebc93f697e38ef0f3aa217)
1dba47a55SKris Buschelman 
2c6db04a5SJed Brown #include <private/pcimpl.h>     /*I "petscpc.h" I*/
33c48a1e8SJed Brown #include <petscdmcomposite.h>   /*I "petscdmcomposite.h" I*/
40971522cSBarry Smith 
5f7c28744SJed Brown const char *const PCFieldSplitSchurPreTypes[] = {"SELF","DIAG","USER","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0};
6f7c28744SJed Brown const char *const PCFieldSplitSchurFactorizationTypes[] = {"DIAG","LOWER","UPPER","FULL","PCFieldSplitSchurFactorizationType","PC_FIELDSPLIT_SCHUR_FACTORIZATION_",0};
7c5d2311dSJed Brown 
8c5d2311dSJed Brown typedef enum {
9c5d2311dSJed Brown   PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG,
10c5d2311dSJed Brown   PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER,
11c5d2311dSJed Brown   PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER,
12c5d2311dSJed Brown   PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL
13c5d2311dSJed Brown } PCFieldSplitSchurFactorizationType;
14084e4875SJed Brown 
150971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
160971522cSBarry Smith struct _PC_FieldSplitLink {
1769a612a9SBarry Smith   KSP               ksp;
180971522cSBarry Smith   Vec               x,y;
19db4c96c1SJed Brown   char              *splitname;
200971522cSBarry Smith   PetscInt          nfields;
210971522cSBarry Smith   PetscInt          *fields;
221b9fc7fcSBarry Smith   VecScatter        sctx;
234aa3045dSJed Brown   IS                is;
2451f519a2SBarry Smith   PC_FieldSplitLink next,previous;
250971522cSBarry Smith };
260971522cSBarry Smith 
270971522cSBarry Smith typedef struct {
2868dd23aaSBarry Smith   PCCompositeType                    type;
29ace3abfcSBarry Smith   PetscBool                          defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */
30ace3abfcSBarry Smith   PetscBool                          splitdefined; /* Flag is set after the splits have been defined, to prevent more splits from being added */
31ace3abfcSBarry Smith   PetscBool                          realdiagonal; /* Flag to use the diagonal blocks of mat preconditioned by pmat, instead of just pmat */
3230ad9308SMatthew Knepley   PetscInt                           bs;           /* Block size for IS and Mat structures */
3330ad9308SMatthew Knepley   PetscInt                           nsplits;      /* Number of field divisions defined */
3479416396SBarry Smith   Vec                                *x,*y,w1,w2;
35519d70e2SJed Brown   Mat                                *mat;         /* The diagonal block for each split */
36519d70e2SJed Brown   Mat                                *pmat;        /* The preconditioning diagonal block for each split */
3730ad9308SMatthew Knepley   Mat                                *Afield;      /* The rows of the matrix associated with each split */
38ace3abfcSBarry Smith   PetscBool                          issetup;
3930ad9308SMatthew Knepley   /* Only used when Schur complement preconditioning is used */
4030ad9308SMatthew Knepley   Mat                                B;            /* The (0,1) block */
4130ad9308SMatthew Knepley   Mat                                C;            /* The (1,0) block */
42a04f6461SBarry Smith   Mat                                schur;        /* The Schur complement S = A11 - A10 A00^{-1} A01 */
43084e4875SJed Brown   Mat                                schur_user;   /* User-provided preconditioning matrix for the Schur complement */
44084e4875SJed Brown   PCFieldSplitSchurPreType           schurpre; /* Determines which preconditioning matrix is used for the Schur complement */
45c5d2311dSJed Brown   PCFieldSplitSchurFactorizationType schurfactorization;
4630ad9308SMatthew Knepley   KSP                                kspschur;     /* The solver for S */
4797bbdb24SBarry Smith   PC_FieldSplitLink                  head;
4863ec74ffSBarry Smith   PetscBool                          reset;         /* indicates PCReset() has been last called on this object, hack */
49c1570756SJed Brown   PetscBool                          suboptionsset; /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */
500971522cSBarry Smith } PC_FieldSplit;
510971522cSBarry Smith 
5216913363SBarry Smith /*
5316913363SBarry Smith     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
5416913363SBarry Smith    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
5516913363SBarry Smith    PC you could change this.
5616913363SBarry Smith */
57084e4875SJed Brown 
58e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it.  This way the
59084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */
60084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac)
61084e4875SJed Brown {
62084e4875SJed Brown   switch (jac->schurpre) {
63084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur;
64084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1];
65084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */
66084e4875SJed Brown     default:
67084e4875SJed Brown       return jac->schur_user ? jac->schur_user : jac->pmat[1];
68084e4875SJed Brown   }
69084e4875SJed Brown }
70084e4875SJed Brown 
71084e4875SJed Brown 
720971522cSBarry Smith #undef __FUNCT__
730971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
740971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
750971522cSBarry Smith {
760971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
770971522cSBarry Smith   PetscErrorCode    ierr;
78ace3abfcSBarry Smith   PetscBool         iascii;
790971522cSBarry Smith   PetscInt          i,j;
805a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
810971522cSBarry Smith 
820971522cSBarry Smith   PetscFunctionBegin;
832692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
840971522cSBarry Smith   if (iascii) {
8551f519a2SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
8669a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
870971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
880971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
891ab39975SBarry Smith       if (ilink->fields) {
900971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
9179416396SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
925a9f2f41SSatish Balay 	for (j=0; j<ilink->nfields; j++) {
9379416396SBarry Smith 	  if (j > 0) {
9479416396SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
9579416396SBarry Smith 	  }
965a9f2f41SSatish Balay 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
970971522cSBarry Smith 	}
980971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
9979416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1001ab39975SBarry Smith       } else {
1011ab39975SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1021ab39975SBarry Smith       }
1035a9f2f41SSatish Balay       ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
1045a9f2f41SSatish Balay       ilink = ilink->next;
1050971522cSBarry Smith     }
1060971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1070971522cSBarry Smith   } else {
10865e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1090971522cSBarry Smith   }
1100971522cSBarry Smith   PetscFunctionReturn(0);
1110971522cSBarry Smith }
1120971522cSBarry Smith 
1130971522cSBarry Smith #undef __FUNCT__
1143b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur"
1153b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
1163b224e63SBarry Smith {
1173b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1183b224e63SBarry Smith   PetscErrorCode    ierr;
119ace3abfcSBarry Smith   PetscBool         iascii;
1203b224e63SBarry Smith   PetscInt          i,j;
1213b224e63SBarry Smith   PC_FieldSplitLink ilink = jac->head;
1223b224e63SBarry Smith   KSP               ksp;
1233b224e63SBarry Smith 
1243b224e63SBarry Smith   PetscFunctionBegin;
1252692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1263b224e63SBarry Smith   if (iascii) {
127c5d2311dSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactorizationTypes[jac->schurfactorization]);CHKERRQ(ierr);
1283b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Split info:\n");CHKERRQ(ierr);
1293b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1303b224e63SBarry Smith     for (i=0; i<jac->nsplits; i++) {
1313b224e63SBarry Smith       if (ilink->fields) {
1323b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
1333b224e63SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
1343b224e63SBarry Smith 	for (j=0; j<ilink->nfields; j++) {
1353b224e63SBarry Smith 	  if (j > 0) {
1363b224e63SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
1373b224e63SBarry Smith 	  }
1383b224e63SBarry Smith 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
1393b224e63SBarry Smith 	}
1403b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1413b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1423b224e63SBarry Smith       } else {
1433b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1443b224e63SBarry Smith       }
1453b224e63SBarry Smith       ilink = ilink->next;
1463b224e63SBarry Smith     }
147435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block \n");CHKERRQ(ierr);
1483b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
14912cae6f2SJed Brown     if (jac->schur) {
15012cae6f2SJed Brown       ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
1513b224e63SBarry Smith       ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
15212cae6f2SJed Brown     } else {
15312cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
15412cae6f2SJed Brown     }
1553b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
156435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr);
1573b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15812cae6f2SJed Brown     if (jac->kspschur) {
1593b224e63SBarry Smith       ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
16012cae6f2SJed Brown     } else {
16112cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
16212cae6f2SJed Brown     }
1633b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1643b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1653b224e63SBarry Smith   } else {
16665e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1673b224e63SBarry Smith   }
1683b224e63SBarry Smith   PetscFunctionReturn(0);
1693b224e63SBarry Smith }
1703b224e63SBarry Smith 
1713b224e63SBarry Smith #undef __FUNCT__
1726c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private"
1736c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */
1746c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc)
1756c924f48SJed Brown {
1766c924f48SJed Brown   PetscErrorCode ierr;
1776c924f48SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1786c924f48SJed Brown   PetscInt       i,nfields,*ifields;
179ace3abfcSBarry Smith   PetscBool      flg;
1806c924f48SJed Brown   char           optionname[128],splitname[8];
1816c924f48SJed Brown 
1826c924f48SJed Brown   PetscFunctionBegin;
1836c924f48SJed Brown   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr);
1846c924f48SJed Brown   for (i=0,flg=PETSC_TRUE; ; i++) {
1856c924f48SJed Brown     ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
1866c924f48SJed Brown     ierr = PetscSNPrintf(optionname,sizeof optionname,"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr);
1876c924f48SJed Brown     nfields = jac->bs;
1886c924f48SJed Brown     ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr);
1896c924f48SJed Brown     if (!flg) break;
1906c924f48SJed Brown     if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
1916c924f48SJed Brown     ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields);CHKERRQ(ierr);
1926c924f48SJed Brown   }
1936c924f48SJed Brown   if (i > 0) {
1946c924f48SJed Brown     /* Makes command-line setting of splits take precedence over setting them in code.
1956c924f48SJed Brown        Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would
1966c924f48SJed Brown        create new splits, which would probably not be what the user wanted. */
1976c924f48SJed Brown     jac->splitdefined = PETSC_TRUE;
1986c924f48SJed Brown   }
1996c924f48SJed Brown   ierr = PetscFree(ifields);CHKERRQ(ierr);
2006c924f48SJed Brown   PetscFunctionReturn(0);
2016c924f48SJed Brown }
2026c924f48SJed Brown 
2036c924f48SJed Brown #undef __FUNCT__
20469a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
20569a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
2060971522cSBarry Smith {
2070971522cSBarry Smith   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
2080971522cSBarry Smith   PetscErrorCode    ierr;
2095a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
2106ce1633cSBarry Smith   PetscBool         flg = PETSC_FALSE,stokes = PETSC_FALSE;
2116c924f48SJed Brown   PetscInt          i;
2120971522cSBarry Smith 
2130971522cSBarry Smith   PetscFunctionBegin;
214d32f9abdSBarry Smith   if (!ilink) {
215d0af7cd3SBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
216d0af7cd3SBarry Smith     if (pc->dm && !stokes) {
2178b8307b2SJed Brown       PetscBool dmcomposite;
2188b8307b2SJed Brown       ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr);
2198b8307b2SJed Brown       if (dmcomposite) {
2208b8307b2SJed Brown         PetscInt nDM;
2218b8307b2SJed Brown         IS       *fields;
2228b8307b2SJed Brown         ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
2238b8307b2SJed Brown         ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr);
2248b8307b2SJed Brown         ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr);
2258b8307b2SJed Brown         for (i=0; i<nDM; i++) {
2268b8307b2SJed Brown           char splitname[8];
2278b8307b2SJed Brown           ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
2288b8307b2SJed Brown           ierr = PCFieldSplitSetIS(pc,splitname,fields[i]);CHKERRQ(ierr);
229fcfd50ebSBarry Smith           ierr = ISDestroy(&fields[i]);CHKERRQ(ierr);
2308b8307b2SJed Brown         }
2318b8307b2SJed Brown         ierr = PetscFree(fields);CHKERRQ(ierr);
2328b8307b2SJed Brown       }
23366ffff09SJed Brown     } else {
234521d7252SBarry Smith       if (jac->bs <= 0) {
235704ba839SBarry Smith         if (pc->pmat) {
236521d7252SBarry Smith           ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
237704ba839SBarry Smith         } else {
238704ba839SBarry Smith           jac->bs = 1;
239704ba839SBarry Smith         }
240521d7252SBarry Smith       }
241d32f9abdSBarry Smith 
242acfcf0e5SJed Brown       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
2436ce1633cSBarry Smith       if (stokes) {
2446ce1633cSBarry Smith         IS       zerodiags,rest;
2456ce1633cSBarry Smith         PetscInt nmin,nmax;
2466ce1633cSBarry Smith 
2476ce1633cSBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
2486ce1633cSBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
2496ce1633cSBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
2506ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
2516ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr);
252fcfd50ebSBarry Smith         ierr = ISDestroy(&zerodiags);CHKERRQ(ierr);
253fcfd50ebSBarry Smith         ierr = ISDestroy(&rest);CHKERRQ(ierr);
2546ce1633cSBarry Smith       } else {
255d32f9abdSBarry Smith         if (!flg) {
256d32f9abdSBarry Smith           /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
257d32f9abdSBarry Smith            then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
2586c924f48SJed Brown           ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
2596c924f48SJed Brown           if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
260d32f9abdSBarry Smith         }
2616c924f48SJed Brown         if (flg || !jac->splitdefined) {
262d32f9abdSBarry Smith           ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
263db4c96c1SJed Brown           for (i=0; i<jac->bs; i++) {
2646c924f48SJed Brown             char splitname[8];
2656c924f48SJed Brown             ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
266db4c96c1SJed Brown             ierr = PCFieldSplitSetFields(pc,splitname,1,&i);CHKERRQ(ierr);
26779416396SBarry Smith           }
26897bbdb24SBarry Smith           jac->defaultsplit = PETSC_TRUE;
269521d7252SBarry Smith         }
27066ffff09SJed Brown       }
2716ce1633cSBarry Smith     }
272edf189efSBarry Smith   } else if (jac->nsplits == 1) {
273edf189efSBarry Smith     if (ilink->is) {
274edf189efSBarry Smith       IS       is2;
275edf189efSBarry Smith       PetscInt nmin,nmax;
276edf189efSBarry Smith 
277edf189efSBarry Smith       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
278edf189efSBarry Smith       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
279db4c96c1SJed Brown       ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr);
280fcfd50ebSBarry Smith       ierr = ISDestroy(&is2);CHKERRQ(ierr);
281e7e72b3dSBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
28263ec74ffSBarry Smith   } else if (jac->reset) {
28363ec74ffSBarry Smith     /* PCReset() has been called on this PC, ilink exists but all IS and Vec data structures in it must be rebuilt
284d0af7cd3SBarry Smith        This is basically the !ilink portion of code above copied from above and the allocation of the ilinks removed
285d0af7cd3SBarry Smith        since they already exist. This should be totally rewritten */
286d0af7cd3SBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
287d0af7cd3SBarry Smith     if (pc->dm && !stokes) {
288d0af7cd3SBarry Smith       PetscBool dmcomposite;
289d0af7cd3SBarry Smith       ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr);
290d0af7cd3SBarry Smith       if (dmcomposite) {
291d0af7cd3SBarry Smith         PetscInt nDM;
292d0af7cd3SBarry Smith         IS       *fields;
293d0af7cd3SBarry Smith         ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
294d0af7cd3SBarry Smith         ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr);
295d0af7cd3SBarry Smith         ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr);
296d0af7cd3SBarry Smith         for (i=0; i<nDM; i++) {
297d0af7cd3SBarry Smith           ilink->is = fields[i];
298d0af7cd3SBarry Smith           ilink     = ilink->next;
299edf189efSBarry Smith         }
300d0af7cd3SBarry Smith         ierr = PetscFree(fields);CHKERRQ(ierr);
301d0af7cd3SBarry Smith       }
302d0af7cd3SBarry Smith     } else {
303d0af7cd3SBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
304d0af7cd3SBarry Smith       if (stokes) {
305d0af7cd3SBarry Smith         IS       zerodiags,rest;
306d0af7cd3SBarry Smith         PetscInt nmin,nmax;
307d0af7cd3SBarry Smith 
308d0af7cd3SBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
309d0af7cd3SBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
310d0af7cd3SBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
31120b26d62SBarry Smith         ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
31220b26d62SBarry Smith         ierr = ISDestroy(&ilink->next->is);CHKERRQ(ierr);
313d0af7cd3SBarry Smith         ilink->is       = rest;
314d0af7cd3SBarry Smith         ilink->next->is = zerodiags;
31520b26d62SBarry Smith       } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used");
316d0af7cd3SBarry Smith     }
317d0af7cd3SBarry Smith   }
318d0af7cd3SBarry Smith 
319e7e72b3dSBarry Smith   if (jac->nsplits < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields");
32069a612a9SBarry Smith   PetscFunctionReturn(0);
32169a612a9SBarry Smith }
32269a612a9SBarry Smith 
32369a612a9SBarry Smith #undef __FUNCT__
32469a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
32569a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
32669a612a9SBarry Smith {
32769a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
32869a612a9SBarry Smith   PetscErrorCode    ierr;
3295a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
330cf502942SBarry Smith   PetscInt          i,nsplit,ccsize;
33169a612a9SBarry Smith   MatStructure      flag = pc->flag;
332ace3abfcSBarry Smith   PetscBool         sorted;
33369a612a9SBarry Smith 
33469a612a9SBarry Smith   PetscFunctionBegin;
33569a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
33697bbdb24SBarry Smith   nsplit = jac->nsplits;
3375a9f2f41SSatish Balay   ilink  = jac->head;
33897bbdb24SBarry Smith 
33997bbdb24SBarry Smith   /* get the matrices for each split */
340704ba839SBarry Smith   if (!jac->issetup) {
3411b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
34297bbdb24SBarry Smith 
343704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
344704ba839SBarry Smith 
345704ba839SBarry Smith     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
34651f519a2SBarry Smith     bs     = jac->bs;
34797bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
348cf502942SBarry Smith     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
3491b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
3501b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
3511b9fc7fcSBarry Smith       if (jac->defaultsplit) {
352704ba839SBarry Smith         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
353704ba839SBarry Smith       } else if (!ilink->is) {
354ccb205f8SBarry Smith         if (ilink->nfields > 1) {
3555a9f2f41SSatish Balay           PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
3565a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
3571b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
3581b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
3591b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
36097bbdb24SBarry Smith             }
36197bbdb24SBarry Smith           }
362d67e408aSBarry Smith           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr);
363ccb205f8SBarry Smith         } else {
364704ba839SBarry Smith           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
365ccb205f8SBarry Smith         }
3663e197d65SBarry Smith       }
367704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
368e32f2f54SBarry Smith       if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
369704ba839SBarry Smith       ilink = ilink->next;
3701b9fc7fcSBarry Smith     }
3711b9fc7fcSBarry Smith   }
3721b9fc7fcSBarry Smith 
373704ba839SBarry Smith   ilink  = jac->head;
37497bbdb24SBarry Smith   if (!jac->pmat) {
375cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
376cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
3774aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
378704ba839SBarry Smith       ilink = ilink->next;
379cf502942SBarry Smith     }
38097bbdb24SBarry Smith   } else {
381cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
3824aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
383704ba839SBarry Smith       ilink = ilink->next;
384cf502942SBarry Smith     }
38597bbdb24SBarry Smith   }
386519d70e2SJed Brown   if (jac->realdiagonal) {
387519d70e2SJed Brown     ilink = jac->head;
388519d70e2SJed Brown     if (!jac->mat) {
389519d70e2SJed Brown       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr);
390519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
391519d70e2SJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
392519d70e2SJed Brown         ilink = ilink->next;
393519d70e2SJed Brown       }
394519d70e2SJed Brown     } else {
395519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
396966e74d7SJed Brown         if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);}
397519d70e2SJed Brown         ilink = ilink->next;
398519d70e2SJed Brown       }
399519d70e2SJed Brown     }
400519d70e2SJed Brown   } else {
401519d70e2SJed Brown     jac->mat = jac->pmat;
402519d70e2SJed Brown   }
40397bbdb24SBarry Smith 
4046c8605c2SJed Brown   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
40568dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
40668dd23aaSBarry Smith     ilink  = jac->head;
40768dd23aaSBarry Smith     if (!jac->Afield) {
40868dd23aaSBarry Smith       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
40968dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
4104aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
41168dd23aaSBarry Smith         ilink = ilink->next;
41268dd23aaSBarry Smith       }
41368dd23aaSBarry Smith     } else {
41468dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
4154aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
41668dd23aaSBarry Smith         ilink = ilink->next;
41768dd23aaSBarry Smith       }
41868dd23aaSBarry Smith     }
41968dd23aaSBarry Smith   }
42068dd23aaSBarry Smith 
4213b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
4223b224e63SBarry Smith     IS       ccis;
4234aa3045dSJed Brown     PetscInt rstart,rend;
424e7e72b3dSBarry Smith     if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
42568dd23aaSBarry Smith 
426e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
427e6cab6aaSJed Brown     ierr  = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
428e6cab6aaSJed Brown 
4293b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
4303b224e63SBarry Smith     if (jac->schur) {
4313b224e63SBarry Smith       ilink = jac->head;
4324aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4334aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
434fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4353b224e63SBarry Smith       ilink = ilink->next;
4364aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4374aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
438fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
439519d70e2SJed Brown       ierr  = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr);
440084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr);
4413b224e63SBarry Smith 
4423b224e63SBarry Smith      } else {
4431cee3971SBarry Smith       KSP ksp;
4446c924f48SJed Brown       char schurprefix[256];
4453b224e63SBarry Smith 
446a04f6461SBarry Smith       /* extract the A01 and A10 matrices */
4473b224e63SBarry Smith       ilink = jac->head;
4484aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4494aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
450fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4513b224e63SBarry Smith       ilink = ilink->next;
4524aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4534aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
454fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4557d50072eSJed Brown       /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */
456519d70e2SJed Brown       ierr  = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr);
457a04f6461SBarry Smith       /* set tabbing and options prefix of KSP inside the MatSchur */
4581cee3971SBarry Smith       ierr  = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
459e69d4d44SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr);
460a04f6461SBarry Smith       ierr  = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",jac->head->splitname);CHKERRQ(ierr);
461a04f6461SBarry Smith       ierr  = KSPSetOptionsPrefix(ksp,schurprefix);CHKERRQ(ierr);
46220b26d62SBarry Smith       /* Need to call this everytime because new matrix is being created */
46320b26d62SBarry Smith       ierr  = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
4643b224e63SBarry Smith 
4653b224e63SBarry Smith       ierr  = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr);
4669005cf84SBarry Smith       ierr  = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr);
4671cee3971SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
468084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
469084e4875SJed Brown       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
470084e4875SJed Brown         PC pc;
471cd5f4a64SJed Brown         ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr);
472084e4875SJed Brown         ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
473084e4875SJed Brown         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
474e69d4d44SBarry Smith       }
4756c924f48SJed Brown       ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
4766c924f48SJed Brown       ierr  = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr);
4773b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
47820b26d62SBarry Smith       /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */
47920b26d62SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
4803b224e63SBarry Smith 
4813b224e63SBarry Smith       ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr);
4823b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr);
4833b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr);
4843b224e63SBarry Smith       ilink = jac->head;
4853b224e63SBarry Smith       ilink->x = jac->x[0]; ilink->y = jac->y[0];
4863b224e63SBarry Smith       ilink = ilink->next;
4873b224e63SBarry Smith       ilink->x = jac->x[1]; ilink->y = jac->y[1];
4883b224e63SBarry Smith     }
4893b224e63SBarry Smith   } else {
49097bbdb24SBarry Smith     /* set up the individual PCs */
49197bbdb24SBarry Smith     i    = 0;
4925a9f2f41SSatish Balay     ilink = jac->head;
4935a9f2f41SSatish Balay     while (ilink) {
494519d70e2SJed Brown       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr);
4953b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
496c1570756SJed Brown       if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);}
4975a9f2f41SSatish Balay       ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
49897bbdb24SBarry Smith       i++;
4995a9f2f41SSatish Balay       ilink = ilink->next;
5000971522cSBarry Smith     }
50197bbdb24SBarry Smith 
50297bbdb24SBarry Smith     /* create work vectors for each split */
5031b9fc7fcSBarry Smith     if (!jac->x) {
50497bbdb24SBarry Smith       ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
5055a9f2f41SSatish Balay       ilink = jac->head;
50697bbdb24SBarry Smith       for (i=0; i<nsplit; i++) {
507906ed7ccSBarry Smith         Vec *vl,*vr;
508906ed7ccSBarry Smith 
509906ed7ccSBarry Smith         ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
510906ed7ccSBarry Smith         ilink->x  = *vr;
511906ed7ccSBarry Smith         ilink->y  = *vl;
512906ed7ccSBarry Smith         ierr      = PetscFree(vr);CHKERRQ(ierr);
513906ed7ccSBarry Smith         ierr      = PetscFree(vl);CHKERRQ(ierr);
5145a9f2f41SSatish Balay         jac->x[i] = ilink->x;
5155a9f2f41SSatish Balay         jac->y[i] = ilink->y;
5165a9f2f41SSatish Balay         ilink     = ilink->next;
51797bbdb24SBarry Smith       }
5183b224e63SBarry Smith     }
5193b224e63SBarry Smith   }
5203b224e63SBarry Smith 
5213b224e63SBarry Smith 
522d0f46423SBarry Smith   if (!jac->head->sctx) {
5233b224e63SBarry Smith     Vec xtmp;
5243b224e63SBarry Smith 
52579416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
5261b9fc7fcSBarry Smith 
5275a9f2f41SSatish Balay     ilink = jac->head;
5281b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
5291b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
530704ba839SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
5315a9f2f41SSatish Balay       ilink = ilink->next;
5321b9fc7fcSBarry Smith     }
5336bf464f9SBarry Smith     ierr = VecDestroy(&xtmp);CHKERRQ(ierr);
5341b9fc7fcSBarry Smith   }
535c1570756SJed Brown   jac->suboptionsset = PETSC_TRUE;
5360971522cSBarry Smith   PetscFunctionReturn(0);
5370971522cSBarry Smith }
5380971522cSBarry Smith 
5395a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
540ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
541ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
5425a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
543ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
544ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
54579416396SBarry Smith 
5460971522cSBarry Smith #undef __FUNCT__
5473b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur"
5483b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
5493b224e63SBarry Smith {
5503b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
5513b224e63SBarry Smith   PetscErrorCode    ierr;
5523b224e63SBarry Smith   KSP               ksp;
5533b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
5543b224e63SBarry Smith 
5553b224e63SBarry Smith   PetscFunctionBegin;
5563b224e63SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
5573b224e63SBarry Smith 
558c5d2311dSJed Brown   switch (jac->schurfactorization) {
559c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG:
560a04f6461SBarry Smith       /* [A00 0; 0 -S], positive definite, suitable for MINRES */
561c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
562c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
563c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
564c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
565c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
566c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
567c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
568c5d2311dSJed Brown       ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr);
569c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
570c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
571c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
572c5d2311dSJed Brown       break;
573c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER:
574a04f6461SBarry Smith       /* [A00 0; A10 S], suitable for left preconditioning */
575c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
576c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
577c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
578c5d2311dSJed Brown       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
579c5d2311dSJed Brown       ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr);
580c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
581c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
582c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
583c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
584c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
585c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
586c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
587c5d2311dSJed Brown       break;
588c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER:
589a04f6461SBarry Smith       /* [A00 A01; 0 S], suitable for right preconditioning */
590c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
591c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
592c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
593c5d2311dSJed Brown       ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
594c5d2311dSJed Brown       ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr);
595c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
596c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
597c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
598c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
599c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
600c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
601c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
602c5d2311dSJed Brown       break;
603c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL:
604a04f6461SBarry 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 */
6053b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6063b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6073b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
6083b224e63SBarry Smith       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
6093b224e63SBarry Smith       ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
6103b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6113b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6123b224e63SBarry Smith 
6133b224e63SBarry Smith       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
6143b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6153b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6163b224e63SBarry Smith 
6173b224e63SBarry Smith       ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
6183b224e63SBarry Smith       ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
6193b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
6203b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6213b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
622c5d2311dSJed Brown   }
6233b224e63SBarry Smith   PetscFunctionReturn(0);
6243b224e63SBarry Smith }
6253b224e63SBarry Smith 
6263b224e63SBarry Smith #undef __FUNCT__
6270971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
6280971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
6290971522cSBarry Smith {
6300971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
6310971522cSBarry Smith   PetscErrorCode    ierr;
6325a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
633af93645bSJed Brown   PetscInt          cnt;
6340971522cSBarry Smith 
6350971522cSBarry Smith   PetscFunctionBegin;
63651f519a2SBarry Smith   CHKMEMQ;
63751f519a2SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
63851f519a2SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
63951f519a2SBarry Smith 
64079416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
6411b9fc7fcSBarry Smith     if (jac->defaultsplit) {
6420971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
6435a9f2f41SSatish Balay       while (ilink) {
6445a9f2f41SSatish Balay         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
6455a9f2f41SSatish Balay         ilink = ilink->next;
6460971522cSBarry Smith       }
6470971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
6481b9fc7fcSBarry Smith     } else {
649efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
6505a9f2f41SSatish Balay       while (ilink) {
6515a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
6525a9f2f41SSatish Balay         ilink = ilink->next;
6531b9fc7fcSBarry Smith       }
6541b9fc7fcSBarry Smith     }
65516913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
65679416396SBarry Smith     if (!jac->w1) {
65779416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
65879416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
65979416396SBarry Smith     }
660efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
6615a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
6623e197d65SBarry Smith     cnt = 1;
6635a9f2f41SSatish Balay     while (ilink->next) {
6645a9f2f41SSatish Balay       ilink = ilink->next;
6653e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
6663e197d65SBarry Smith       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
6673e197d65SBarry Smith       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
6683e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6693e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6703e197d65SBarry Smith       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
6713e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6723e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6733e197d65SBarry Smith     }
67451f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
67511755939SBarry Smith       cnt -= 2;
67651f519a2SBarry Smith       while (ilink->previous) {
67751f519a2SBarry Smith         ilink = ilink->previous;
67811755939SBarry Smith         /* compute the residual only over the part of the vector needed */
67911755939SBarry Smith         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
68011755939SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
68111755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
68211755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
68311755939SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
68411755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
68511755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
68651f519a2SBarry Smith       }
68711755939SBarry Smith     }
68865e19b50SBarry Smith   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
68951f519a2SBarry Smith   CHKMEMQ;
6900971522cSBarry Smith   PetscFunctionReturn(0);
6910971522cSBarry Smith }
6920971522cSBarry Smith 
693421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
694ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
695ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
696421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
697ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
698ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
699421e10b8SBarry Smith 
700421e10b8SBarry Smith #undef __FUNCT__
7018c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit"
702421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
703421e10b8SBarry Smith {
704421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
705421e10b8SBarry Smith   PetscErrorCode    ierr;
706421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
707421e10b8SBarry Smith 
708421e10b8SBarry Smith   PetscFunctionBegin;
709421e10b8SBarry Smith   CHKMEMQ;
710421e10b8SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
711421e10b8SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
712421e10b8SBarry Smith 
713421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
714421e10b8SBarry Smith     if (jac->defaultsplit) {
715421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
716421e10b8SBarry Smith       while (ilink) {
717421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
718421e10b8SBarry Smith 	ilink = ilink->next;
719421e10b8SBarry Smith       }
720421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
721421e10b8SBarry Smith     } else {
722421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
723421e10b8SBarry Smith       while (ilink) {
724421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
725421e10b8SBarry Smith 	ilink = ilink->next;
726421e10b8SBarry Smith       }
727421e10b8SBarry Smith     }
728421e10b8SBarry Smith   } else {
729421e10b8SBarry Smith     if (!jac->w1) {
730421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
731421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
732421e10b8SBarry Smith     }
733421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
734421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
735421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
736421e10b8SBarry Smith       while (ilink->next) {
737421e10b8SBarry Smith         ilink = ilink->next;
7389989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
739421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
740421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
741421e10b8SBarry Smith       }
742421e10b8SBarry Smith       while (ilink->previous) {
743421e10b8SBarry Smith         ilink = ilink->previous;
7449989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
745421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
746421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
747421e10b8SBarry Smith       }
748421e10b8SBarry Smith     } else {
749421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
750421e10b8SBarry Smith 	ilink = ilink->next;
751421e10b8SBarry Smith       }
752421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
753421e10b8SBarry Smith       while (ilink->previous) {
754421e10b8SBarry Smith 	ilink = ilink->previous;
7559989ab13SBarry Smith 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
756421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
757421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
758421e10b8SBarry Smith       }
759421e10b8SBarry Smith     }
760421e10b8SBarry Smith   }
761421e10b8SBarry Smith   CHKMEMQ;
762421e10b8SBarry Smith   PetscFunctionReturn(0);
763421e10b8SBarry Smith }
764421e10b8SBarry Smith 
7650971522cSBarry Smith #undef __FUNCT__
766574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit"
767574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc)
7680971522cSBarry Smith {
7690971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
7700971522cSBarry Smith   PetscErrorCode    ierr;
7715a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
7720971522cSBarry Smith 
7730971522cSBarry Smith   PetscFunctionBegin;
7745a9f2f41SSatish Balay   while (ilink) {
775574deadeSBarry Smith     ierr = KSPReset(ilink->ksp);CHKERRQ(ierr);
776fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->x);CHKERRQ(ierr);
777fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->y);CHKERRQ(ierr);
778fcfd50ebSBarry Smith     ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr);
779fcfd50ebSBarry Smith     ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
7805a9f2f41SSatish Balay     next = ilink->next;
7815a9f2f41SSatish Balay     ilink = next;
7820971522cSBarry Smith   }
78305b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
784574deadeSBarry Smith   if (jac->mat && jac->mat != jac->pmat) {
785574deadeSBarry Smith     ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);
786574deadeSBarry Smith   } else if (jac->mat) {
787574deadeSBarry Smith     jac->mat = PETSC_NULL;
788574deadeSBarry Smith   }
78997bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
79068dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
7916bf464f9SBarry Smith   ierr = VecDestroy(&jac->w1);CHKERRQ(ierr);
7926bf464f9SBarry Smith   ierr = VecDestroy(&jac->w2);CHKERRQ(ierr);
7936bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur);CHKERRQ(ierr);
7946bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
7956bf464f9SBarry Smith   ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr);
7966bf464f9SBarry Smith   ierr = MatDestroy(&jac->B);CHKERRQ(ierr);
7976bf464f9SBarry Smith   ierr = MatDestroy(&jac->C);CHKERRQ(ierr);
79863ec74ffSBarry Smith   jac->reset = PETSC_TRUE;
799574deadeSBarry Smith   PetscFunctionReturn(0);
800574deadeSBarry Smith }
801574deadeSBarry Smith 
802574deadeSBarry Smith #undef __FUNCT__
803574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
804574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
805574deadeSBarry Smith {
806574deadeSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
807574deadeSBarry Smith   PetscErrorCode    ierr;
808574deadeSBarry Smith   PC_FieldSplitLink ilink = jac->head,next;
809574deadeSBarry Smith 
810574deadeSBarry Smith   PetscFunctionBegin;
811574deadeSBarry Smith   ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr);
812574deadeSBarry Smith   while (ilink) {
8136bf464f9SBarry Smith     ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr);
814574deadeSBarry Smith     next = ilink->next;
815574deadeSBarry Smith     ierr = PetscFree(ilink->splitname);CHKERRQ(ierr);
816574deadeSBarry Smith     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
817574deadeSBarry Smith     ierr = PetscFree(ilink);CHKERRQ(ierr);
818574deadeSBarry Smith     ilink = next;
819574deadeSBarry Smith   }
820574deadeSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
821c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
8220971522cSBarry Smith   PetscFunctionReturn(0);
8230971522cSBarry Smith }
8240971522cSBarry Smith 
8250971522cSBarry Smith #undef __FUNCT__
8260971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
8270971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
8280971522cSBarry Smith {
8291b9fc7fcSBarry Smith   PetscErrorCode  ierr;
8306c924f48SJed Brown   PetscInt        bs;
831bc59fbc5SBarry Smith   PetscBool       flg,stokes = PETSC_FALSE;
8329dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
8333b224e63SBarry Smith   PCCompositeType ctype;
8341b9fc7fcSBarry Smith 
8350971522cSBarry Smith   PetscFunctionBegin;
8361b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
837acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr);
83851f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
83951f519a2SBarry Smith   if (flg) {
84051f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
84151f519a2SBarry Smith   }
842704ba839SBarry Smith 
843435f959eSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
844c0adfefeSBarry Smith   if (stokes) {
845c0adfefeSBarry Smith     ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
846c0adfefeSBarry Smith     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
847c0adfefeSBarry Smith   }
848c0adfefeSBarry Smith 
8493b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
8503b224e63SBarry Smith   if (flg) {
8513b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
8523b224e63SBarry Smith   }
853d32f9abdSBarry Smith 
854c30613efSMatthew Knepley   /* Only setup fields once */
855c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
856d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
857d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
8586c924f48SJed Brown     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
8596c924f48SJed Brown     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
860d32f9abdSBarry Smith   }
861c5d2311dSJed Brown   if (jac->type == PC_COMPOSITE_SCHUR) {
862c5d2311dSJed Brown     ierr = PetscOptionsEnum("-pc_fieldsplit_schur_factorization_type","Factorization to use","None",PCFieldSplitSchurFactorizationTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,PETSC_NULL);CHKERRQ(ierr);
863084e4875SJed 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);
864c5d2311dSJed Brown   }
8651b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
8660971522cSBarry Smith   PetscFunctionReturn(0);
8670971522cSBarry Smith }
8680971522cSBarry Smith 
8690971522cSBarry Smith /*------------------------------------------------------------------------------------*/
8700971522cSBarry Smith 
8710971522cSBarry Smith EXTERN_C_BEGIN
8720971522cSBarry Smith #undef __FUNCT__
8730971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
8747087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields)
8750971522cSBarry Smith {
87697bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
8770971522cSBarry Smith   PetscErrorCode    ierr;
8785a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
87969a612a9SBarry Smith   char              prefix[128];
88051f519a2SBarry Smith   PetscInt          i;
8810971522cSBarry Smith 
8820971522cSBarry Smith   PetscFunctionBegin;
8836c924f48SJed Brown   if (jac->splitdefined) {
8846c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
8856c924f48SJed Brown     PetscFunctionReturn(0);
8866c924f48SJed Brown   }
88751f519a2SBarry Smith   for (i=0; i<n; i++) {
888e32f2f54SBarry 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);
889e32f2f54SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
89051f519a2SBarry Smith   }
891704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
892a04f6461SBarry Smith   if (splitname) {
893db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
894a04f6461SBarry Smith   } else {
895a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
896a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
897a04f6461SBarry Smith   }
898704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
8995a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
9005a9f2f41SSatish Balay   ilink->nfields = n;
9015a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
9027adad957SLisandro Dalcin   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
9031cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
9045a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
9059005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
90669a612a9SBarry Smith 
907a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
9085a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
9090971522cSBarry Smith 
9100971522cSBarry Smith   if (!next) {
9115a9f2f41SSatish Balay     jac->head       = ilink;
91251f519a2SBarry Smith     ilink->previous = PETSC_NULL;
9130971522cSBarry Smith   } else {
9140971522cSBarry Smith     while (next->next) {
9150971522cSBarry Smith       next = next->next;
9160971522cSBarry Smith     }
9175a9f2f41SSatish Balay     next->next      = ilink;
91851f519a2SBarry Smith     ilink->previous = next;
9190971522cSBarry Smith   }
9200971522cSBarry Smith   jac->nsplits++;
9210971522cSBarry Smith   PetscFunctionReturn(0);
9220971522cSBarry Smith }
9230971522cSBarry Smith EXTERN_C_END
9240971522cSBarry Smith 
925e69d4d44SBarry Smith EXTERN_C_BEGIN
926e69d4d44SBarry Smith #undef __FUNCT__
927e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
9287087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
929e69d4d44SBarry Smith {
930e69d4d44SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
931e69d4d44SBarry Smith   PetscErrorCode ierr;
932e69d4d44SBarry Smith 
933e69d4d44SBarry Smith   PetscFunctionBegin;
934519d70e2SJed Brown   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
935e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
936e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
93713e0d083SBarry Smith   if (n) *n = jac->nsplits;
938e69d4d44SBarry Smith   PetscFunctionReturn(0);
939e69d4d44SBarry Smith }
940e69d4d44SBarry Smith EXTERN_C_END
9410971522cSBarry Smith 
9420971522cSBarry Smith EXTERN_C_BEGIN
9430971522cSBarry Smith #undef __FUNCT__
94469a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
9457087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
9460971522cSBarry Smith {
9470971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
9480971522cSBarry Smith   PetscErrorCode    ierr;
9490971522cSBarry Smith   PetscInt          cnt = 0;
9505a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
9510971522cSBarry Smith 
9520971522cSBarry Smith   PetscFunctionBegin;
9535d480477SMatthew G Knepley   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
9545a9f2f41SSatish Balay   while (ilink) {
9555a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
9565a9f2f41SSatish Balay     ilink = ilink->next;
9570971522cSBarry Smith   }
9585d480477SMatthew 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);
95913e0d083SBarry Smith   if (n) *n = jac->nsplits;
9600971522cSBarry Smith   PetscFunctionReturn(0);
9610971522cSBarry Smith }
9620971522cSBarry Smith EXTERN_C_END
9630971522cSBarry Smith 
964704ba839SBarry Smith EXTERN_C_BEGIN
965704ba839SBarry Smith #undef __FUNCT__
966704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
9677087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
968704ba839SBarry Smith {
969704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
970704ba839SBarry Smith   PetscErrorCode    ierr;
971704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
972704ba839SBarry Smith   char              prefix[128];
973704ba839SBarry Smith 
974704ba839SBarry Smith   PetscFunctionBegin;
9756c924f48SJed Brown   if (jac->splitdefined) {
9766c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
9776c924f48SJed Brown     PetscFunctionReturn(0);
9786c924f48SJed Brown   }
97916913363SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
980a04f6461SBarry Smith   if (splitname) {
981db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
982a04f6461SBarry Smith   } else {
983a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
984a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
985a04f6461SBarry Smith   }
9861ab39975SBarry Smith   ilink->is      = is;
9871ab39975SBarry Smith   ierr           = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
988704ba839SBarry Smith   ilink->next    = PETSC_NULL;
989704ba839SBarry Smith   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
9901cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
991704ba839SBarry Smith   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
9929005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
993704ba839SBarry Smith 
994a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
995704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
996704ba839SBarry Smith 
997704ba839SBarry Smith   if (!next) {
998704ba839SBarry Smith     jac->head       = ilink;
999704ba839SBarry Smith     ilink->previous = PETSC_NULL;
1000704ba839SBarry Smith   } else {
1001704ba839SBarry Smith     while (next->next) {
1002704ba839SBarry Smith       next = next->next;
1003704ba839SBarry Smith     }
1004704ba839SBarry Smith     next->next      = ilink;
1005704ba839SBarry Smith     ilink->previous = next;
1006704ba839SBarry Smith   }
1007704ba839SBarry Smith   jac->nsplits++;
1008704ba839SBarry Smith 
1009704ba839SBarry Smith   PetscFunctionReturn(0);
1010704ba839SBarry Smith }
1011704ba839SBarry Smith EXTERN_C_END
1012704ba839SBarry Smith 
10130971522cSBarry Smith #undef __FUNCT__
10140971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
10150971522cSBarry Smith /*@
10160971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
10170971522cSBarry Smith 
1018ad4df100SBarry Smith     Logically Collective on PC
10190971522cSBarry Smith 
10200971522cSBarry Smith     Input Parameters:
10210971522cSBarry Smith +   pc  - the preconditioner context
1022a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
10230971522cSBarry Smith .   n - the number of fields in this split
1024db4c96c1SJed Brown -   fields - the fields in this split
10250971522cSBarry Smith 
10260971522cSBarry Smith     Level: intermediate
10270971522cSBarry Smith 
1028d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
1029d32f9abdSBarry Smith 
1030d32f9abdSBarry Smith      The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block
1031d32f9abdSBarry 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
1032d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1033d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
1034d32f9abdSBarry Smith 
1035db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
1036db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1037db4c96c1SJed Brown 
1038d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
10390971522cSBarry Smith 
10400971522cSBarry Smith @*/
10417087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields)
10420971522cSBarry Smith {
10434ac538c5SBarry Smith   PetscErrorCode ierr;
10440971522cSBarry Smith 
10450971522cSBarry Smith   PetscFunctionBegin;
10460700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1047db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1048db4c96c1SJed Brown   if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1049db4c96c1SJed Brown   PetscValidIntPointer(fields,3);
10504ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *),(pc,splitname,n,fields));CHKERRQ(ierr);
10510971522cSBarry Smith   PetscFunctionReturn(0);
10520971522cSBarry Smith }
10530971522cSBarry Smith 
10540971522cSBarry Smith #undef __FUNCT__
1055704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
1056704ba839SBarry Smith /*@
1057704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
1058704ba839SBarry Smith 
1059ad4df100SBarry Smith     Logically Collective on PC
1060704ba839SBarry Smith 
1061704ba839SBarry Smith     Input Parameters:
1062704ba839SBarry Smith +   pc  - the preconditioner context
1063a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
1064db4c96c1SJed Brown -   is - the index set that defines the vector elements in this field
1065704ba839SBarry Smith 
1066d32f9abdSBarry Smith 
1067a6ffb8dbSJed Brown     Notes:
1068a6ffb8dbSJed Brown     Use PCFieldSplitSetFields(), for fields defined by strided types.
1069a6ffb8dbSJed Brown 
1070db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
1071db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1072d32f9abdSBarry Smith 
1073704ba839SBarry Smith     Level: intermediate
1074704ba839SBarry Smith 
1075704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
1076704ba839SBarry Smith 
1077704ba839SBarry Smith @*/
10787087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1079704ba839SBarry Smith {
10804ac538c5SBarry Smith   PetscErrorCode ierr;
1081704ba839SBarry Smith 
1082704ba839SBarry Smith   PetscFunctionBegin;
10830700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1084db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1085db4c96c1SJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,3);
10864ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr);
1087704ba839SBarry Smith   PetscFunctionReturn(0);
1088704ba839SBarry Smith }
1089704ba839SBarry Smith 
1090704ba839SBarry Smith #undef __FUNCT__
109157a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS"
109257a9adfeSMatthew G Knepley /*@
109357a9adfeSMatthew G Knepley     PCFieldSplitGetIS - Retrieves the elements for a field as an IS
109457a9adfeSMatthew G Knepley 
109557a9adfeSMatthew G Knepley     Logically Collective on PC
109657a9adfeSMatthew G Knepley 
109757a9adfeSMatthew G Knepley     Input Parameters:
109857a9adfeSMatthew G Knepley +   pc  - the preconditioner context
109957a9adfeSMatthew G Knepley -   splitname - name of this split
110057a9adfeSMatthew G Knepley 
110157a9adfeSMatthew G Knepley     Output Parameter:
110257a9adfeSMatthew G Knepley -   is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found
110357a9adfeSMatthew G Knepley 
110457a9adfeSMatthew G Knepley     Level: intermediate
110557a9adfeSMatthew G Knepley 
110657a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS()
110757a9adfeSMatthew G Knepley 
110857a9adfeSMatthew G Knepley @*/
110957a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
111057a9adfeSMatthew G Knepley {
111157a9adfeSMatthew G Knepley   PetscErrorCode ierr;
111257a9adfeSMatthew G Knepley 
111357a9adfeSMatthew G Knepley   PetscFunctionBegin;
111457a9adfeSMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
111557a9adfeSMatthew G Knepley   PetscValidCharPointer(splitname,2);
111657a9adfeSMatthew G Knepley   PetscValidPointer(is,3);
111757a9adfeSMatthew G Knepley   {
111857a9adfeSMatthew G Knepley     PC_FieldSplit    *jac   = (PC_FieldSplit *) pc->data;
111957a9adfeSMatthew G Knepley     PC_FieldSplitLink ilink = jac->head;
112057a9adfeSMatthew G Knepley     PetscBool         found;
112157a9adfeSMatthew G Knepley 
112257a9adfeSMatthew G Knepley     *is = PETSC_NULL;
112357a9adfeSMatthew G Knepley     while(ilink) {
112457a9adfeSMatthew G Knepley       ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr);
112557a9adfeSMatthew G Knepley       if (found) {
112657a9adfeSMatthew G Knepley         *is = ilink->is;
112757a9adfeSMatthew G Knepley         break;
112857a9adfeSMatthew G Knepley       }
112957a9adfeSMatthew G Knepley       ilink = ilink->next;
113057a9adfeSMatthew G Knepley     }
113157a9adfeSMatthew G Knepley   }
113257a9adfeSMatthew G Knepley   PetscFunctionReturn(0);
113357a9adfeSMatthew G Knepley }
113457a9adfeSMatthew G Knepley 
113557a9adfeSMatthew G Knepley #undef __FUNCT__
113651f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
113751f519a2SBarry Smith /*@
113851f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
113951f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
114051f519a2SBarry Smith 
1141ad4df100SBarry Smith     Logically Collective on PC
114251f519a2SBarry Smith 
114351f519a2SBarry Smith     Input Parameters:
114451f519a2SBarry Smith +   pc  - the preconditioner context
114551f519a2SBarry Smith -   bs - the block size
114651f519a2SBarry Smith 
114751f519a2SBarry Smith     Level: intermediate
114851f519a2SBarry Smith 
114951f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
115051f519a2SBarry Smith 
115151f519a2SBarry Smith @*/
11527087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
115351f519a2SBarry Smith {
11544ac538c5SBarry Smith   PetscErrorCode ierr;
115551f519a2SBarry Smith 
115651f519a2SBarry Smith   PetscFunctionBegin;
11570700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1158c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,bs,2);
11594ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr);
116051f519a2SBarry Smith   PetscFunctionReturn(0);
116151f519a2SBarry Smith }
116251f519a2SBarry Smith 
116351f519a2SBarry Smith #undef __FUNCT__
116469a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
11650971522cSBarry Smith /*@C
116669a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
11670971522cSBarry Smith 
116869a612a9SBarry Smith    Collective on KSP
11690971522cSBarry Smith 
11700971522cSBarry Smith    Input Parameter:
11710971522cSBarry Smith .  pc - the preconditioner context
11720971522cSBarry Smith 
11730971522cSBarry Smith    Output Parameters:
117413e0d083SBarry Smith +  n - the number of splits
117569a612a9SBarry Smith -  pc - the array of KSP contexts
11760971522cSBarry Smith 
11770971522cSBarry Smith    Note:
1178d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
1179d32f9abdSBarry Smith    (not the KSP just the array that contains them).
11800971522cSBarry Smith 
118169a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
11820971522cSBarry Smith 
11830971522cSBarry Smith    Level: advanced
11840971522cSBarry Smith 
11850971522cSBarry Smith .seealso: PCFIELDSPLIT
11860971522cSBarry Smith @*/
11877087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
11880971522cSBarry Smith {
11894ac538c5SBarry Smith   PetscErrorCode ierr;
11900971522cSBarry Smith 
11910971522cSBarry Smith   PetscFunctionBegin;
11920700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
119313e0d083SBarry Smith   if (n) PetscValidIntPointer(n,2);
11944ac538c5SBarry Smith   ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
11950971522cSBarry Smith   PetscFunctionReturn(0);
11960971522cSBarry Smith }
11970971522cSBarry Smith 
1198e69d4d44SBarry Smith #undef __FUNCT__
1199e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition"
1200e69d4d44SBarry Smith /*@
1201e69d4d44SBarry Smith     PCFieldSplitSchurPrecondition -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
1202a04f6461SBarry Smith       A11 matrix. Otherwise no preconditioner is used.
1203e69d4d44SBarry Smith 
1204e69d4d44SBarry Smith     Collective on PC
1205e69d4d44SBarry Smith 
1206e69d4d44SBarry Smith     Input Parameters:
1207e69d4d44SBarry Smith +   pc  - the preconditioner context
1208fd1303e9SJungho Lee .   ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default
1209084e4875SJed Brown -   userpre - matrix to use for preconditioning, or PETSC_NULL
1210084e4875SJed Brown 
1211e69d4d44SBarry Smith     Options Database:
1212084e4875SJed Brown .     -pc_fieldsplit_schur_precondition <self,user,diag> default is diag
1213e69d4d44SBarry Smith 
1214fd1303e9SJungho Lee     Notes:
1215fd1303e9SJungho Lee $    If ptype is
1216fd1303e9SJungho Lee $        user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument
1217fd1303e9SJungho Lee $             to this function).
1218fd1303e9SJungho Lee $        diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original
1219fd1303e9SJungho Lee $             matrix associated with the Schur complement (i.e. A11)
1220fd1303e9SJungho Lee $        self the preconditioner for the Schur complement is generated from the Schur complement matrix itself:
1221fd1303e9SJungho Lee $             The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC
1222fd1303e9SJungho Lee $             preconditioner
1223fd1303e9SJungho Lee 
1224fd1303e9SJungho Lee      When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense
1225fd1303e9SJungho Lee     with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
1226fd1303e9SJungho Lee     -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement.
1227fd1303e9SJungho Lee 
1228fd1303e9SJungho Lee     Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in
1229fd1303e9SJungho Lee      the name since it sets a proceedure to use.
1230fd1303e9SJungho Lee 
1231e69d4d44SBarry Smith     Level: intermediate
1232e69d4d44SBarry Smith 
1233fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC
1234e69d4d44SBarry Smith 
1235e69d4d44SBarry Smith @*/
12367087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1237e69d4d44SBarry Smith {
12384ac538c5SBarry Smith   PetscErrorCode ierr;
1239e69d4d44SBarry Smith 
1240e69d4d44SBarry Smith   PetscFunctionBegin;
12410700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12424ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr);
1243e69d4d44SBarry Smith   PetscFunctionReturn(0);
1244e69d4d44SBarry Smith }
1245e69d4d44SBarry Smith 
1246e69d4d44SBarry Smith EXTERN_C_BEGIN
1247e69d4d44SBarry Smith #undef __FUNCT__
1248e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
12497087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1250e69d4d44SBarry Smith {
1251e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1252084e4875SJed Brown   PetscErrorCode  ierr;
1253e69d4d44SBarry Smith 
1254e69d4d44SBarry Smith   PetscFunctionBegin;
1255084e4875SJed Brown   jac->schurpre = ptype;
1256084e4875SJed Brown   if (pre) {
12576bf464f9SBarry Smith     ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
1258084e4875SJed Brown     jac->schur_user = pre;
1259084e4875SJed Brown     ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1260084e4875SJed Brown   }
1261e69d4d44SBarry Smith   PetscFunctionReturn(0);
1262e69d4d44SBarry Smith }
1263e69d4d44SBarry Smith EXTERN_C_END
1264e69d4d44SBarry Smith 
126530ad9308SMatthew Knepley #undef __FUNCT__
126630ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
126730ad9308SMatthew Knepley /*@C
12688c03b21aSDmitry Karpeev    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
126930ad9308SMatthew Knepley 
127030ad9308SMatthew Knepley    Collective on KSP
127130ad9308SMatthew Knepley 
127230ad9308SMatthew Knepley    Input Parameter:
127330ad9308SMatthew Knepley .  pc - the preconditioner context
127430ad9308SMatthew Knepley 
127530ad9308SMatthew Knepley    Output Parameters:
1276a04f6461SBarry Smith +  A00 - the (0,0) block
1277a04f6461SBarry Smith .  A01 - the (0,1) block
1278a04f6461SBarry Smith .  A10 - the (1,0) block
1279a04f6461SBarry Smith -  A11 - the (1,1) block
128030ad9308SMatthew Knepley 
128130ad9308SMatthew Knepley    Level: advanced
128230ad9308SMatthew Knepley 
128330ad9308SMatthew Knepley .seealso: PCFIELDSPLIT
128430ad9308SMatthew Knepley @*/
1285a04f6461SBarry Smith PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
128630ad9308SMatthew Knepley {
128730ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit *) pc->data;
128830ad9308SMatthew Knepley 
128930ad9308SMatthew Knepley   PetscFunctionBegin;
12900700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1291c1235816SBarry Smith   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
1292a04f6461SBarry Smith   if (A00) *A00 = jac->pmat[0];
1293a04f6461SBarry Smith   if (A01) *A01 = jac->B;
1294a04f6461SBarry Smith   if (A10) *A10 = jac->C;
1295a04f6461SBarry Smith   if (A11) *A11 = jac->pmat[1];
129630ad9308SMatthew Knepley   PetscFunctionReturn(0);
129730ad9308SMatthew Knepley }
129830ad9308SMatthew Knepley 
129979416396SBarry Smith EXTERN_C_BEGIN
130079416396SBarry Smith #undef __FUNCT__
130179416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
13027087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
130379416396SBarry Smith {
130479416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1305e69d4d44SBarry Smith   PetscErrorCode ierr;
130679416396SBarry Smith 
130779416396SBarry Smith   PetscFunctionBegin;
130879416396SBarry Smith   jac->type = type;
13093b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
13103b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
13113b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
1312e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1313e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1314e69d4d44SBarry Smith 
13153b224e63SBarry Smith   } else {
13163b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
13173b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
1318e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
13199e7d6b0aSBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr);
13203b224e63SBarry Smith   }
132179416396SBarry Smith   PetscFunctionReturn(0);
132279416396SBarry Smith }
132379416396SBarry Smith EXTERN_C_END
132479416396SBarry Smith 
132551f519a2SBarry Smith EXTERN_C_BEGIN
132651f519a2SBarry Smith #undef __FUNCT__
132751f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
13287087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
132951f519a2SBarry Smith {
133051f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
133151f519a2SBarry Smith 
133251f519a2SBarry Smith   PetscFunctionBegin;
133365e19b50SBarry Smith   if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
133465e19b50SBarry 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);
133551f519a2SBarry Smith   jac->bs = bs;
133651f519a2SBarry Smith   PetscFunctionReturn(0);
133751f519a2SBarry Smith }
133851f519a2SBarry Smith EXTERN_C_END
133951f519a2SBarry Smith 
134079416396SBarry Smith #undef __FUNCT__
134179416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
1342bc08b0f1SBarry Smith /*@
134379416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
134479416396SBarry Smith 
134579416396SBarry Smith    Collective on PC
134679416396SBarry Smith 
134779416396SBarry Smith    Input Parameter:
134879416396SBarry Smith .  pc - the preconditioner context
134981540f2fSBarry Smith .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
135079416396SBarry Smith 
135179416396SBarry Smith    Options Database Key:
1352a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
135379416396SBarry Smith 
135479416396SBarry Smith    Level: Developer
135579416396SBarry Smith 
135679416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
135779416396SBarry Smith 
135879416396SBarry Smith .seealso: PCCompositeSetType()
135979416396SBarry Smith 
136079416396SBarry Smith @*/
13617087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
136279416396SBarry Smith {
13634ac538c5SBarry Smith   PetscErrorCode ierr;
136479416396SBarry Smith 
136579416396SBarry Smith   PetscFunctionBegin;
13660700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13674ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr);
136879416396SBarry Smith   PetscFunctionReturn(0);
136979416396SBarry Smith }
137079416396SBarry Smith 
13710971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
13720971522cSBarry Smith /*MC
1373a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
1374a04f6461SBarry Smith                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.
13750971522cSBarry Smith 
1376edf189efSBarry Smith      To set options on the solvers for each block append -fieldsplit_ to all the PC
1377edf189efSBarry Smith         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
13780971522cSBarry Smith 
1379a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
138069a612a9SBarry Smith          and set the options directly on the resulting KSP object
13810971522cSBarry Smith 
13820971522cSBarry Smith    Level: intermediate
13830971522cSBarry Smith 
138479416396SBarry Smith    Options Database Keys:
138581540f2fSBarry Smith +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
138681540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
138781540f2fSBarry Smith                               been supplied explicitly by -pc_fieldsplit_%d_fields
138881540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
138981540f2fSBarry Smith .   -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative>
1390e69d4d44SBarry Smith .   -pc_fieldsplit_schur_precondition <true,false> default is true
1391435f959eSBarry Smith .   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver
139279416396SBarry Smith 
1393edf189efSBarry Smith -    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
13943b224e63SBarry Smith      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
13953b224e63SBarry Smith 
1396d32f9abdSBarry Smith    Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1397d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
1398d32f9abdSBarry Smith 
1399d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
1400d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1401d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
1402d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
1403d32f9abdSBarry Smith 
1404a04f6461SBarry Smith      For the Schur complement preconditioner if J = ( A00 A01 )
1405a04f6461SBarry Smith                                                     ( A10 A11 )
1406e69d4d44SBarry Smith      the preconditioner is
1407a04f6461SBarry Smith               (I   -B ksp(A00)) ( inv(A00)   0  (I         0  )
1408a04f6461SBarry Smith               (0    I       ) (   0    ksp(S) ) (-A10 ksp(A00) I  )
1409a04f6461SBarry Smith      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of
1410a04f6461SBarry 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).
1411a04f6461SBarry Smith      For PCFieldSplitGetKSP() when field number is
1412edf189efSBarry Smith      0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
1413a04f6461SBarry Smith      A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
14147e8cb189SBarry Smith      option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc
1415e69d4d44SBarry Smith 
1416edf189efSBarry Smith      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
1417edf189efSBarry Smith      is used automatically for a second block.
1418edf189efSBarry Smith 
1419*ff218e97SBarry Smith      The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1.
1420*ff218e97SBarry Smith      Generally it should be used with the AIJ format.
1421*ff218e97SBarry Smith 
1422*ff218e97SBarry Smith      The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see,
1423*ff218e97SBarry Smith      for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT
1424*ff218e97SBarry Smith      inside a smoother resulting in "Distributive Smoothers".
14250716a85fSBarry Smith 
1426a541d17aSBarry Smith    Concepts: physics based preconditioners, block preconditioners
14270971522cSBarry Smith 
14287e8cb189SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
1429e69d4d44SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
14300971522cSBarry Smith M*/
14310971522cSBarry Smith 
14320971522cSBarry Smith EXTERN_C_BEGIN
14330971522cSBarry Smith #undef __FUNCT__
14340971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
14357087cfbeSBarry Smith PetscErrorCode  PCCreate_FieldSplit(PC pc)
14360971522cSBarry Smith {
14370971522cSBarry Smith   PetscErrorCode ierr;
14380971522cSBarry Smith   PC_FieldSplit  *jac;
14390971522cSBarry Smith 
14400971522cSBarry Smith   PetscFunctionBegin;
144138f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
14420971522cSBarry Smith   jac->bs        = -1;
14430971522cSBarry Smith   jac->nsplits   = 0;
14443e197d65SBarry Smith   jac->type      = PC_COMPOSITE_MULTIPLICATIVE;
1445e6cab6aaSJed Brown   jac->schurpre  = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
1446c5d2311dSJed Brown   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL;
144751f519a2SBarry Smith 
14480971522cSBarry Smith   pc->data     = (void*)jac;
14490971522cSBarry Smith 
14500971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
1451421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
14520971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
1453574deadeSBarry Smith   pc->ops->reset             = PCReset_FieldSplit;
14540971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
14550971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
14560971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
14570971522cSBarry Smith   pc->ops->applyrichardson   = 0;
14580971522cSBarry Smith 
145969a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
146069a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
14610971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
14620971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1463704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
1464704ba839SBarry Smith                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
146579416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
146679416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
146751f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
146851f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
14690971522cSBarry Smith   PetscFunctionReturn(0);
14700971522cSBarry Smith }
14710971522cSBarry Smith EXTERN_C_END
14720971522cSBarry Smith 
14730971522cSBarry Smith 
1474a541d17aSBarry Smith 
1475