xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 4d343eea47f3e457863a63f23fb9efcfef73ea2c)
1dba47a55SKris Buschelman 
2c6db04a5SJed Brown #include <private/pcimpl.h>     /*I "petscpc.h" I*/
33c48a1e8SJed Brown #include <petscdmcomposite.h>   /*I "petscdmcomposite.h" I*/
4d2822287SMatthew G Knepley #include <petscdmcomplex.h>
50971522cSBarry Smith 
6f7c28744SJed Brown const char *const PCFieldSplitSchurPreTypes[] = {"SELF","DIAG","USER","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0};
7f7c28744SJed Brown const char *const PCFieldSplitSchurFactorizationTypes[] = {"DIAG","LOWER","UPPER","FULL","PCFieldSplitSchurFactorizationType","PC_FIELDSPLIT_SCHUR_FACTORIZATION_",0};
8c5d2311dSJed Brown 
9c5d2311dSJed Brown typedef enum {
10c5d2311dSJed Brown   PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG,
11c5d2311dSJed Brown   PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER,
12c5d2311dSJed Brown   PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER,
13c5d2311dSJed Brown   PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL
14c5d2311dSJed Brown } PCFieldSplitSchurFactorizationType;
15084e4875SJed Brown 
160971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
170971522cSBarry Smith struct _PC_FieldSplitLink {
1869a612a9SBarry Smith   KSP               ksp;
190971522cSBarry Smith   Vec               x,y;
20db4c96c1SJed Brown   char              *splitname;
210971522cSBarry Smith   PetscInt          nfields;
220971522cSBarry Smith   PetscInt          *fields;
231b9fc7fcSBarry Smith   VecScatter        sctx;
244aa3045dSJed Brown   IS                is;
2551f519a2SBarry Smith   PC_FieldSplitLink next,previous;
260971522cSBarry Smith };
270971522cSBarry Smith 
280971522cSBarry Smith typedef struct {
2968dd23aaSBarry Smith   PCCompositeType                    type;
30ace3abfcSBarry Smith   PetscBool                          defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */
31ace3abfcSBarry Smith   PetscBool                          splitdefined; /* Flag is set after the splits have been defined, to prevent more splits from being added */
32ace3abfcSBarry Smith   PetscBool                          realdiagonal; /* Flag to use the diagonal blocks of mat preconditioned by pmat, instead of just pmat */
3330ad9308SMatthew Knepley   PetscInt                           bs;           /* Block size for IS and Mat structures */
3430ad9308SMatthew Knepley   PetscInt                           nsplits;      /* Number of field divisions defined */
3579416396SBarry Smith   Vec                                *x,*y,w1,w2;
36519d70e2SJed Brown   Mat                                *mat;         /* The diagonal block for each split */
37519d70e2SJed Brown   Mat                                *pmat;        /* The preconditioning diagonal block for each split */
3830ad9308SMatthew Knepley   Mat                                *Afield;      /* The rows of the matrix associated with each split */
39ace3abfcSBarry Smith   PetscBool                          issetup;
4030ad9308SMatthew Knepley   /* Only used when Schur complement preconditioning is used */
4130ad9308SMatthew Knepley   Mat                                B;            /* The (0,1) block */
4230ad9308SMatthew Knepley   Mat                                C;            /* The (1,0) block */
43a04f6461SBarry Smith   Mat                                schur;        /* The Schur complement S = A11 - A10 A00^{-1} A01 */
44084e4875SJed Brown   Mat                                schur_user;   /* User-provided preconditioning matrix for the Schur complement */
45084e4875SJed Brown   PCFieldSplitSchurPreType           schurpre; /* Determines which preconditioning matrix is used for the Schur complement */
46c5d2311dSJed Brown   PCFieldSplitSchurFactorizationType schurfactorization;
4730ad9308SMatthew Knepley   KSP                                kspschur;     /* The solver for S */
4897bbdb24SBarry Smith   PC_FieldSplitLink                  head;
4963ec74ffSBarry Smith   PetscBool                          reset;         /* indicates PCReset() has been last called on this object, hack */
50c1570756SJed Brown   PetscBool                          suboptionsset; /* Indicates that the KSPSetFromOptions() has been called on the sub-KSPs */
510971522cSBarry Smith } PC_FieldSplit;
520971522cSBarry Smith 
5316913363SBarry Smith /*
5416913363SBarry Smith     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
5516913363SBarry Smith    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
5616913363SBarry Smith    PC you could change this.
5716913363SBarry Smith */
58084e4875SJed Brown 
59e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it.  This way the
60084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */
61084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac)
62084e4875SJed Brown {
63084e4875SJed Brown   switch (jac->schurpre) {
64084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur;
65084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1];
66084e4875SJed Brown     case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */
67084e4875SJed Brown     default:
68084e4875SJed Brown       return jac->schur_user ? jac->schur_user : jac->pmat[1];
69084e4875SJed Brown   }
70084e4875SJed Brown }
71084e4875SJed Brown 
72084e4875SJed Brown 
730971522cSBarry Smith #undef __FUNCT__
740971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
750971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
760971522cSBarry Smith {
770971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
780971522cSBarry Smith   PetscErrorCode    ierr;
79ace3abfcSBarry Smith   PetscBool         iascii;
800971522cSBarry Smith   PetscInt          i,j;
815a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
820971522cSBarry Smith 
830971522cSBarry Smith   PetscFunctionBegin;
842692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
850971522cSBarry Smith   if (iascii) {
8651f519a2SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
8769a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
880971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
890971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
901ab39975SBarry Smith       if (ilink->fields) {
910971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
9279416396SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
935a9f2f41SSatish Balay 	for (j=0; j<ilink->nfields; j++) {
9479416396SBarry Smith 	  if (j > 0) {
9579416396SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
9679416396SBarry Smith 	  }
975a9f2f41SSatish Balay 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
980971522cSBarry Smith 	}
990971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
10079416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1011ab39975SBarry Smith       } else {
1021ab39975SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1031ab39975SBarry Smith       }
1045a9f2f41SSatish Balay       ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
1055a9f2f41SSatish Balay       ilink = ilink->next;
1060971522cSBarry Smith     }
1070971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1080971522cSBarry Smith   } else {
10965e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1100971522cSBarry Smith   }
1110971522cSBarry Smith   PetscFunctionReturn(0);
1120971522cSBarry Smith }
1130971522cSBarry Smith 
1140971522cSBarry Smith #undef __FUNCT__
1153b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur"
1163b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
1173b224e63SBarry Smith {
1183b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
1193b224e63SBarry Smith   PetscErrorCode    ierr;
120ace3abfcSBarry Smith   PetscBool         iascii;
1213b224e63SBarry Smith   PetscInt          i,j;
1223b224e63SBarry Smith   PC_FieldSplitLink ilink = jac->head;
1233b224e63SBarry Smith   KSP               ksp;
1243b224e63SBarry Smith 
1253b224e63SBarry Smith   PetscFunctionBegin;
1262692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1273b224e63SBarry Smith   if (iascii) {
128c5d2311dSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, blocksize = %D, factorization %s\n",jac->bs,PCFieldSplitSchurFactorizationTypes[jac->schurfactorization]);CHKERRQ(ierr);
1293b224e63SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Split info:\n");CHKERRQ(ierr);
1303b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1313b224e63SBarry Smith     for (i=0; i<jac->nsplits; i++) {
1323b224e63SBarry Smith       if (ilink->fields) {
1333b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
1343b224e63SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
1353b224e63SBarry Smith 	for (j=0; j<ilink->nfields; j++) {
1363b224e63SBarry Smith 	  if (j > 0) {
1373b224e63SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
1383b224e63SBarry Smith 	  }
1393b224e63SBarry Smith 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
1403b224e63SBarry Smith 	}
1413b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1423b224e63SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1433b224e63SBarry Smith       } else {
1443b224e63SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
1453b224e63SBarry Smith       }
1463b224e63SBarry Smith       ilink = ilink->next;
1473b224e63SBarry Smith     }
148435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A00 block \n");CHKERRQ(ierr);
1493b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15012cae6f2SJed Brown     if (jac->schur) {
15112cae6f2SJed Brown       ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
1523b224e63SBarry Smith       ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
15312cae6f2SJed Brown     } else {
15412cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
15512cae6f2SJed Brown     }
1563b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
157435f959eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = A11 - A10 inv(A00) A01 \n");CHKERRQ(ierr);
1583b224e63SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15912cae6f2SJed Brown     if (jac->kspschur) {
1603b224e63SBarry Smith       ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
16112cae6f2SJed Brown     } else {
16212cae6f2SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
16312cae6f2SJed Brown     }
1643b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1653b224e63SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1663b224e63SBarry Smith   } else {
16765e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
1683b224e63SBarry Smith   }
1693b224e63SBarry Smith   PetscFunctionReturn(0);
1703b224e63SBarry Smith }
1713b224e63SBarry Smith 
1723b224e63SBarry Smith #undef __FUNCT__
1736c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private"
1746c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */
1756c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc)
1766c924f48SJed Brown {
1776c924f48SJed Brown   PetscErrorCode ierr;
1786c924f48SJed Brown   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1796c924f48SJed Brown   PetscInt       i,nfields,*ifields;
180ace3abfcSBarry Smith   PetscBool      flg;
1816c924f48SJed Brown   char           optionname[128],splitname[8];
1826c924f48SJed Brown 
1836c924f48SJed Brown   PetscFunctionBegin;
1846c924f48SJed Brown   ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr);
1856c924f48SJed Brown   for (i=0,flg=PETSC_TRUE; ; i++) {
1866c924f48SJed Brown     ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
1876c924f48SJed Brown     ierr = PetscSNPrintf(optionname,sizeof optionname,"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr);
1886c924f48SJed Brown     nfields = jac->bs;
1896c924f48SJed Brown     ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr);
1906c924f48SJed Brown     if (!flg) break;
191ea6813d4SMatthew G Knepley     if (!nfields) SETERRQ(((PetscObject) pc)->comm,PETSC_ERR_USER,"Cannot list zero fields");
1926c924f48SJed Brown     ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields);CHKERRQ(ierr);
1936c924f48SJed Brown   }
1946c924f48SJed Brown   if (i > 0) {
1956c924f48SJed Brown     /* Makes command-line setting of splits take precedence over setting them in code.
1966c924f48SJed Brown        Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would
1976c924f48SJed Brown        create new splits, which would probably not be what the user wanted. */
1986c924f48SJed Brown     jac->splitdefined = PETSC_TRUE;
1996c924f48SJed Brown   }
2006c924f48SJed Brown   ierr = PetscFree(ifields);CHKERRQ(ierr);
2016c924f48SJed Brown   PetscFunctionReturn(0);
2026c924f48SJed Brown }
2036c924f48SJed Brown 
2046c924f48SJed Brown #undef __FUNCT__
20569a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
20669a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
2070971522cSBarry Smith {
2080971522cSBarry Smith   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
2090971522cSBarry Smith   PetscErrorCode    ierr;
2105a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
2116ce1633cSBarry Smith   PetscBool         flg = PETSC_FALSE,stokes = PETSC_FALSE;
2126c924f48SJed Brown   PetscInt          i;
2130971522cSBarry Smith 
2140971522cSBarry Smith   PetscFunctionBegin;
215d32f9abdSBarry Smith   if (!ilink) {
216d0af7cd3SBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
217d0af7cd3SBarry Smith     if (pc->dm && !stokes) {
218*4d343eeaSMatthew G Knepley       PetscInt     numFields, f;
219*4d343eeaSMatthew G Knepley       const char **fieldNames;
2208b8307b2SJed Brown       IS          *fields;
221*4d343eeaSMatthew G Knepley       PetscBool    dmcomposite;
222*4d343eeaSMatthew G Knepley 
223*4d343eeaSMatthew G Knepley       ierr = DMCreateFieldIS(pc->dm, &numFields, &fieldNames, &fields);CHKERRQ(ierr);
224*4d343eeaSMatthew G Knepley       for(f = 0; f < numFields; ++f) {
225*4d343eeaSMatthew G Knepley         ierr = PCFieldSplitSetIS(pc, fieldNames[f], fields[f]);CHKERRQ(ierr);
226*4d343eeaSMatthew G Knepley         ierr = PetscFree(fieldNames[f]);CHKERRQ(ierr);
227*4d343eeaSMatthew G Knepley         ierr = ISDestroy(&fields[f]);CHKERRQ(ierr);
228*4d343eeaSMatthew G Knepley       }
229*4d343eeaSMatthew G Knepley       ierr = PetscFree(fieldNames);CHKERRQ(ierr);
230*4d343eeaSMatthew G Knepley       ierr = PetscFree(fields);CHKERRQ(ierr);
231*4d343eeaSMatthew G Knepley 
232*4d343eeaSMatthew G Knepley       ierr = PetscTypeCompare((PetscObject) pc->dm, DMCOMPOSITE, &dmcomposite);CHKERRQ(ierr);
233*4d343eeaSMatthew G Knepley       if (dmcomposite) {
2342fa5ba8aSJed Brown         DM      *dms;
235*4d343eeaSMatthew G Knepley         PetscInt nDM;
236*4d343eeaSMatthew G Knepley 
2378b8307b2SJed Brown         ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
2388b8307b2SJed Brown         ierr = DMCompositeGetNumberDM(pc->dm, &nDM);CHKERRQ(ierr);
2392fa5ba8aSJed Brown         ierr = PetscMalloc(nDM*sizeof(DM), &dms);CHKERRQ(ierr);
2402fa5ba8aSJed Brown         ierr = DMCompositeGetEntriesArray(pc->dm, dms);CHKERRQ(ierr);
2412fa5ba8aSJed Brown         for (ilink=jac->head,i=0; ilink; ilink=ilink->next,i++) {
2422fa5ba8aSJed Brown           ierr = KSPSetDM(ilink->ksp, dms[i]);CHKERRQ(ierr);
2432fa5ba8aSJed Brown           ierr = KSPSetDMActive(ilink->ksp, PETSC_FALSE);CHKERRQ(ierr);
2442fa5ba8aSJed Brown         }
2452fa5ba8aSJed Brown         ierr = PetscFree(dms);CHKERRQ(ierr);
2468b8307b2SJed Brown       }
24766ffff09SJed Brown     } else {
248521d7252SBarry Smith       if (jac->bs <= 0) {
249704ba839SBarry Smith         if (pc->pmat) {
250521d7252SBarry Smith           ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
251704ba839SBarry Smith         } else {
252704ba839SBarry Smith           jac->bs = 1;
253704ba839SBarry Smith         }
254521d7252SBarry Smith       }
255d32f9abdSBarry Smith 
256acfcf0e5SJed Brown       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
2576ce1633cSBarry Smith       if (stokes) {
2586ce1633cSBarry Smith         IS       zerodiags,rest;
2596ce1633cSBarry Smith         PetscInt nmin,nmax;
2606ce1633cSBarry Smith 
2616ce1633cSBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
2626ce1633cSBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
2636ce1633cSBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
2646ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"0",rest);CHKERRQ(ierr);
2656ce1633cSBarry Smith         ierr = PCFieldSplitSetIS(pc,"1",zerodiags);CHKERRQ(ierr);
266fcfd50ebSBarry Smith         ierr = ISDestroy(&zerodiags);CHKERRQ(ierr);
267fcfd50ebSBarry Smith         ierr = ISDestroy(&rest);CHKERRQ(ierr);
2686ce1633cSBarry Smith       } else {
269d32f9abdSBarry Smith         if (!flg) {
270d32f9abdSBarry Smith           /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
271d32f9abdSBarry Smith            then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
2726c924f48SJed Brown           ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
2736c924f48SJed Brown           if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
274d32f9abdSBarry Smith         }
2756c924f48SJed Brown         if (flg || !jac->splitdefined) {
276d32f9abdSBarry Smith           ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
277db4c96c1SJed Brown           for (i=0; i<jac->bs; i++) {
2786c924f48SJed Brown             char splitname[8];
2796c924f48SJed Brown             ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr);
280db4c96c1SJed Brown             ierr = PCFieldSplitSetFields(pc,splitname,1,&i);CHKERRQ(ierr);
28179416396SBarry Smith           }
28297bbdb24SBarry Smith           jac->defaultsplit = PETSC_TRUE;
283521d7252SBarry Smith         }
28466ffff09SJed Brown       }
2856ce1633cSBarry Smith     }
286edf189efSBarry Smith   } else if (jac->nsplits == 1) {
287edf189efSBarry Smith     if (ilink->is) {
288edf189efSBarry Smith       IS       is2;
289edf189efSBarry Smith       PetscInt nmin,nmax;
290edf189efSBarry Smith 
291edf189efSBarry Smith       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
292edf189efSBarry Smith       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
293db4c96c1SJed Brown       ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr);
294fcfd50ebSBarry Smith       ierr = ISDestroy(&is2);CHKERRQ(ierr);
295ea6813d4SMatthew G Knepley     } else SETERRQ(((PetscObject) pc)->comm,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
29663ec74ffSBarry Smith   } else if (jac->reset) {
29763ec74ffSBarry Smith     /* PCReset() has been called on this PC, ilink exists but all IS and Vec data structures in it must be rebuilt
298d0af7cd3SBarry Smith        This is basically the !ilink portion of code above copied from above and the allocation of the ilinks removed
299d0af7cd3SBarry Smith        since they already exist. This should be totally rewritten */
300d0af7cd3SBarry Smith     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
301d0af7cd3SBarry Smith     if (pc->dm && !stokes) {
302d0af7cd3SBarry Smith       PetscBool dmcomposite;
303d0af7cd3SBarry Smith       ierr = PetscTypeCompare((PetscObject)pc->dm,DMCOMPOSITE,&dmcomposite);CHKERRQ(ierr);
304d0af7cd3SBarry Smith       if (dmcomposite) {
305d0af7cd3SBarry Smith         PetscInt nDM;
306d0af7cd3SBarry Smith         IS       *fields;
3072fa5ba8aSJed Brown         DM       *dms;
308d0af7cd3SBarry Smith         ierr = PetscInfo(pc,"Setting up physics based fieldsplit preconditioner using the embedded DM\n");CHKERRQ(ierr);
309d0af7cd3SBarry Smith         ierr = DMCompositeGetNumberDM(pc->dm,&nDM);CHKERRQ(ierr);
310d0af7cd3SBarry Smith         ierr = DMCompositeGetGlobalISs(pc->dm,&fields);CHKERRQ(ierr);
3112fa5ba8aSJed Brown         ierr = PetscMalloc(nDM*sizeof(DM),&dms);CHKERRQ(ierr);
3122fa5ba8aSJed Brown         ierr = DMCompositeGetEntriesArray(pc->dm,dms);CHKERRQ(ierr);
313d0af7cd3SBarry Smith         for (i=0; i<nDM; i++) {
3142fa5ba8aSJed Brown           ierr = KSPSetDM(ilink->ksp,dms[i]);CHKERRQ(ierr);
3152fa5ba8aSJed Brown           ierr = KSPSetDMActive(ilink->ksp,PETSC_FALSE);CHKERRQ(ierr);
316d0af7cd3SBarry Smith           ilink->is = fields[i];
317d0af7cd3SBarry Smith           ilink     = ilink->next;
318edf189efSBarry Smith         }
319d0af7cd3SBarry Smith         ierr = PetscFree(fields);CHKERRQ(ierr);
3202fa5ba8aSJed Brown         ierr = PetscFree(dms);CHKERRQ(ierr);
321d0af7cd3SBarry Smith       }
322d0af7cd3SBarry Smith     } else {
323d0af7cd3SBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
324d0af7cd3SBarry Smith       if (stokes) {
325d0af7cd3SBarry Smith         IS       zerodiags,rest;
326d0af7cd3SBarry Smith         PetscInt nmin,nmax;
327d0af7cd3SBarry Smith 
328d0af7cd3SBarry Smith         ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
329d0af7cd3SBarry Smith         ierr = MatFindZeroDiagonals(pc->mat,&zerodiags);CHKERRQ(ierr);
330d0af7cd3SBarry Smith         ierr = ISComplement(zerodiags,nmin,nmax,&rest);CHKERRQ(ierr);
33120b26d62SBarry Smith         ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
33220b26d62SBarry Smith         ierr = ISDestroy(&ilink->next->is);CHKERRQ(ierr);
333d0af7cd3SBarry Smith         ilink->is       = rest;
334d0af7cd3SBarry Smith         ilink->next->is = zerodiags;
33520b26d62SBarry Smith       } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cases not yet handled when PCReset() was used");
336d0af7cd3SBarry Smith     }
337d0af7cd3SBarry Smith   }
338d0af7cd3SBarry Smith 
339ea6813d4SMatthew G Knepley   if (jac->nsplits < 2) SETERRQ1(((PetscObject) pc)->comm,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields, not %d", jac->nsplits);
34069a612a9SBarry Smith   PetscFunctionReturn(0);
34169a612a9SBarry Smith }
34269a612a9SBarry Smith 
34369a612a9SBarry Smith #undef __FUNCT__
34469a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
34569a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
34669a612a9SBarry Smith {
34769a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
34869a612a9SBarry Smith   PetscErrorCode    ierr;
3495a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
350cf502942SBarry Smith   PetscInt          i,nsplit,ccsize;
35169a612a9SBarry Smith   MatStructure      flag = pc->flag;
352ace3abfcSBarry Smith   PetscBool         sorted;
35369a612a9SBarry Smith 
35469a612a9SBarry Smith   PetscFunctionBegin;
35569a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
35697bbdb24SBarry Smith   nsplit = jac->nsplits;
3575a9f2f41SSatish Balay   ilink  = jac->head;
35897bbdb24SBarry Smith 
35997bbdb24SBarry Smith   /* get the matrices for each split */
360704ba839SBarry Smith   if (!jac->issetup) {
3611b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
36297bbdb24SBarry Smith 
363704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
364704ba839SBarry Smith 
365704ba839SBarry Smith     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
36651f519a2SBarry Smith     bs     = jac->bs;
36797bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
368cf502942SBarry Smith     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
3691b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
3701b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
3711b9fc7fcSBarry Smith       if (jac->defaultsplit) {
372704ba839SBarry Smith         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
373704ba839SBarry Smith       } else if (!ilink->is) {
374ccb205f8SBarry Smith         if (ilink->nfields > 1) {
3755a9f2f41SSatish Balay           PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
3765a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
3771b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
3781b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
3791b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
38097bbdb24SBarry Smith             }
38197bbdb24SBarry Smith           }
382d67e408aSBarry Smith           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,PETSC_OWN_POINTER,&ilink->is);CHKERRQ(ierr);
383ccb205f8SBarry Smith         } else {
384704ba839SBarry Smith           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
385ccb205f8SBarry Smith         }
3863e197d65SBarry Smith       }
387704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
388e32f2f54SBarry Smith       if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
389704ba839SBarry Smith       ilink = ilink->next;
3901b9fc7fcSBarry Smith     }
3911b9fc7fcSBarry Smith   }
3921b9fc7fcSBarry Smith 
393704ba839SBarry Smith   ilink  = jac->head;
39497bbdb24SBarry Smith   if (!jac->pmat) {
395cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
396cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
3974aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
398704ba839SBarry Smith       ilink = ilink->next;
399cf502942SBarry Smith     }
40097bbdb24SBarry Smith   } else {
401cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
4024aa3045dSJed Brown       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
403704ba839SBarry Smith       ilink = ilink->next;
404cf502942SBarry Smith     }
40597bbdb24SBarry Smith   }
406519d70e2SJed Brown   if (jac->realdiagonal) {
407519d70e2SJed Brown     ilink = jac->head;
408519d70e2SJed Brown     if (!jac->mat) {
409519d70e2SJed Brown       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr);
410519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
411519d70e2SJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
412519d70e2SJed Brown         ilink = ilink->next;
413519d70e2SJed Brown       }
414519d70e2SJed Brown     } else {
415519d70e2SJed Brown       for (i=0; i<nsplit; i++) {
416966e74d7SJed Brown         if (jac->mat[i]) {ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);}
417519d70e2SJed Brown         ilink = ilink->next;
418519d70e2SJed Brown       }
419519d70e2SJed Brown     }
420519d70e2SJed Brown   } else {
421519d70e2SJed Brown     jac->mat = jac->pmat;
422519d70e2SJed Brown   }
42397bbdb24SBarry Smith 
4246c8605c2SJed Brown   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
42568dd23aaSBarry Smith     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
42668dd23aaSBarry Smith     ilink  = jac->head;
42768dd23aaSBarry Smith     if (!jac->Afield) {
42868dd23aaSBarry Smith       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
42968dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
4304aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
43168dd23aaSBarry Smith         ilink = ilink->next;
43268dd23aaSBarry Smith       }
43368dd23aaSBarry Smith     } else {
43468dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
4354aa3045dSJed Brown         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
43668dd23aaSBarry Smith         ilink = ilink->next;
43768dd23aaSBarry Smith       }
43868dd23aaSBarry Smith     }
43968dd23aaSBarry Smith   }
44068dd23aaSBarry Smith 
4413b224e63SBarry Smith   if (jac->type == PC_COMPOSITE_SCHUR) {
4423b224e63SBarry Smith     IS       ccis;
4434aa3045dSJed Brown     PetscInt rstart,rend;
444e7e72b3dSBarry Smith     if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
44568dd23aaSBarry Smith 
446e6cab6aaSJed Brown     /* When extracting off-diagonal submatrices, we take complements from this range */
447e6cab6aaSJed Brown     ierr  = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
448e6cab6aaSJed Brown 
4493b224e63SBarry Smith     /* need to handle case when one is resetting up the preconditioner */
4503b224e63SBarry Smith     if (jac->schur) {
4513b224e63SBarry Smith       ilink = jac->head;
4524aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4534aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
454fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4553b224e63SBarry Smith       ilink = ilink->next;
4564aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4574aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
458fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
459519d70e2SJed Brown       ierr  = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr);
460084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr);
4613b224e63SBarry Smith 
4623b224e63SBarry Smith      } else {
4631cee3971SBarry Smith       KSP ksp;
4646c924f48SJed Brown       char schurprefix[256];
4653b224e63SBarry Smith 
466a04f6461SBarry Smith       /* extract the A01 and A10 matrices */
4673b224e63SBarry Smith       ilink = jac->head;
4684aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4694aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
470fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4713b224e63SBarry Smith       ilink = ilink->next;
4724aa3045dSJed Brown       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
4734aa3045dSJed Brown       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
474fcfd50ebSBarry Smith       ierr  = ISDestroy(&ccis);CHKERRQ(ierr);
4757d50072eSJed Brown       /* Use mat[0] (diagonal block of the real matrix) preconditioned by pmat[0] */
476519d70e2SJed Brown       ierr  = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr);
477a04f6461SBarry Smith       /* set tabbing and options prefix of KSP inside the MatSchur */
4781cee3971SBarry Smith       ierr  = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
479e69d4d44SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr);
480a04f6461SBarry Smith       ierr  = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",jac->head->splitname);CHKERRQ(ierr);
481a04f6461SBarry Smith       ierr  = KSPSetOptionsPrefix(ksp,schurprefix);CHKERRQ(ierr);
48220b26d62SBarry Smith       /* Need to call this everytime because new matrix is being created */
48320b26d62SBarry Smith       ierr  = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
4847ae8954aSSatish Balay       ierr  = MatSetUp(jac->schur);CHKERRQ(ierr);
4853b224e63SBarry Smith 
4863b224e63SBarry Smith       ierr  = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr);
4879005cf84SBarry Smith       ierr  = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->kspschur);CHKERRQ(ierr);
4881cee3971SBarry Smith       ierr  = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
489084e4875SJed Brown       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
490084e4875SJed Brown       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
491084e4875SJed Brown         PC pc;
492cd5f4a64SJed Brown         ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr);
493084e4875SJed Brown         ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
494084e4875SJed Brown         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
495e69d4d44SBarry Smith       }
4966c924f48SJed Brown       ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
4976c924f48SJed Brown       ierr  = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr);
4983b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
49920b26d62SBarry Smith       /* need to call this every time, since the jac->kspschur is freshly created, otherwise its options never get set */
50020b26d62SBarry Smith       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
5013b224e63SBarry Smith 
5023b224e63SBarry Smith       ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr);
5033b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr);
5043b224e63SBarry Smith       ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr);
5053b224e63SBarry Smith       ilink = jac->head;
5063b224e63SBarry Smith       ilink->x = jac->x[0]; ilink->y = jac->y[0];
5073b224e63SBarry Smith       ilink = ilink->next;
5083b224e63SBarry Smith       ilink->x = jac->x[1]; ilink->y = jac->y[1];
5093b224e63SBarry Smith     }
5103b224e63SBarry Smith   } else {
51197bbdb24SBarry Smith     /* set up the individual PCs */
51297bbdb24SBarry Smith     i    = 0;
5135a9f2f41SSatish Balay     ilink = jac->head;
5145a9f2f41SSatish Balay     while (ilink) {
515519d70e2SJed Brown       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr);
5163b224e63SBarry Smith       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
517c1570756SJed Brown       if (!jac->suboptionsset) {ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);}
5185a9f2f41SSatish Balay       ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
51997bbdb24SBarry Smith       i++;
5205a9f2f41SSatish Balay       ilink = ilink->next;
5210971522cSBarry Smith     }
52297bbdb24SBarry Smith 
52397bbdb24SBarry Smith     /* create work vectors for each split */
5241b9fc7fcSBarry Smith     if (!jac->x) {
52597bbdb24SBarry Smith       ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
5265a9f2f41SSatish Balay       ilink = jac->head;
52797bbdb24SBarry Smith       for (i=0; i<nsplit; i++) {
528906ed7ccSBarry Smith         Vec *vl,*vr;
529906ed7ccSBarry Smith 
530906ed7ccSBarry Smith         ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
531906ed7ccSBarry Smith         ilink->x  = *vr;
532906ed7ccSBarry Smith         ilink->y  = *vl;
533906ed7ccSBarry Smith         ierr      = PetscFree(vr);CHKERRQ(ierr);
534906ed7ccSBarry Smith         ierr      = PetscFree(vl);CHKERRQ(ierr);
5355a9f2f41SSatish Balay         jac->x[i] = ilink->x;
5365a9f2f41SSatish Balay         jac->y[i] = ilink->y;
5375a9f2f41SSatish Balay         ilink     = ilink->next;
53897bbdb24SBarry Smith       }
5393b224e63SBarry Smith     }
5403b224e63SBarry Smith   }
5413b224e63SBarry Smith 
5423b224e63SBarry Smith 
543d0f46423SBarry Smith   if (!jac->head->sctx) {
5443b224e63SBarry Smith     Vec xtmp;
5453b224e63SBarry Smith 
54679416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
5471b9fc7fcSBarry Smith 
5485a9f2f41SSatish Balay     ilink = jac->head;
5491b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
5501b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
551704ba839SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
5525a9f2f41SSatish Balay       ilink = ilink->next;
5531b9fc7fcSBarry Smith     }
5546bf464f9SBarry Smith     ierr = VecDestroy(&xtmp);CHKERRQ(ierr);
5551b9fc7fcSBarry Smith   }
556c1570756SJed Brown   jac->suboptionsset = PETSC_TRUE;
5570971522cSBarry Smith   PetscFunctionReturn(0);
5580971522cSBarry Smith }
5590971522cSBarry Smith 
5605a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
561ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
562ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
5635a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
564ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
565ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
56679416396SBarry Smith 
5670971522cSBarry Smith #undef __FUNCT__
5683b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur"
5693b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
5703b224e63SBarry Smith {
5713b224e63SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
5723b224e63SBarry Smith   PetscErrorCode    ierr;
5733b224e63SBarry Smith   KSP               ksp;
5743b224e63SBarry Smith   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
5753b224e63SBarry Smith 
5763b224e63SBarry Smith   PetscFunctionBegin;
5773b224e63SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
5783b224e63SBarry Smith 
579c5d2311dSJed Brown   switch (jac->schurfactorization) {
580c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_DIAG:
581a04f6461SBarry Smith       /* [A00 0; 0 -S], positive definite, suitable for MINRES */
582c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
583c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
584c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
585c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
586c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
587c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
588c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
589c5d2311dSJed Brown       ierr = VecScale(ilinkD->y,-1.);CHKERRQ(ierr);
590c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
591c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
592c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
593c5d2311dSJed Brown       break;
594c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_LOWER:
595a04f6461SBarry Smith       /* [A00 0; A10 S], suitable for left preconditioning */
596c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
597c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
598c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
599c5d2311dSJed Brown       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
600c5d2311dSJed Brown       ierr = VecScale(ilinkD->x,-1.);CHKERRQ(ierr);
601c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
602c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
603c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
604c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
605c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
606c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
607c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
608c5d2311dSJed Brown       break;
609c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_UPPER:
610a04f6461SBarry Smith       /* [A00 A01; 0 S], suitable for right preconditioning */
611c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
612c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
613c5d2311dSJed Brown       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
614c5d2311dSJed Brown       ierr = MatMult(jac->B,ilinkD->y,ilinkA->x);CHKERRQ(ierr);
615c5d2311dSJed Brown       ierr = VecScale(ilinkA->x,-1.);CHKERRQ(ierr);
616c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
617c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
618c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
619c5d2311dSJed Brown       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
620c5d2311dSJed Brown       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
621c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
622c5d2311dSJed Brown       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
623c5d2311dSJed Brown       break;
624c5d2311dSJed Brown     case PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL:
625a04f6461SBarry 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 */
6263b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6273b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6283b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
6293b224e63SBarry Smith       ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
6303b224e63SBarry Smith       ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
6313b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6323b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6333b224e63SBarry Smith 
6343b224e63SBarry Smith       ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
6353b224e63SBarry Smith       ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6363b224e63SBarry Smith       ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6373b224e63SBarry Smith 
6383b224e63SBarry Smith       ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
6393b224e63SBarry Smith       ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
6403b224e63SBarry Smith       ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
6413b224e63SBarry Smith       ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6423b224e63SBarry Smith       ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
643c5d2311dSJed Brown   }
6443b224e63SBarry Smith   PetscFunctionReturn(0);
6453b224e63SBarry Smith }
6463b224e63SBarry Smith 
6473b224e63SBarry Smith #undef __FUNCT__
6480971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
6490971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
6500971522cSBarry Smith {
6510971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
6520971522cSBarry Smith   PetscErrorCode    ierr;
6535a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
654af93645bSJed Brown   PetscInt          cnt;
6550971522cSBarry Smith 
6560971522cSBarry Smith   PetscFunctionBegin;
65751f519a2SBarry Smith   CHKMEMQ;
65851f519a2SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
65951f519a2SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
66051f519a2SBarry Smith 
66179416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
6621b9fc7fcSBarry Smith     if (jac->defaultsplit) {
6630971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
6645a9f2f41SSatish Balay       while (ilink) {
6655a9f2f41SSatish Balay         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
6665a9f2f41SSatish Balay         ilink = ilink->next;
6670971522cSBarry Smith       }
6680971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
6691b9fc7fcSBarry Smith     } else {
670efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
6715a9f2f41SSatish Balay       while (ilink) {
6725a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
6735a9f2f41SSatish Balay         ilink = ilink->next;
6741b9fc7fcSBarry Smith       }
6751b9fc7fcSBarry Smith     }
67616913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
67779416396SBarry Smith     if (!jac->w1) {
67879416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
67979416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
68079416396SBarry Smith     }
681efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
6825a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
6833e197d65SBarry Smith     cnt = 1;
6845a9f2f41SSatish Balay     while (ilink->next) {
6855a9f2f41SSatish Balay       ilink = ilink->next;
6863e197d65SBarry Smith       /* compute the residual only over the part of the vector needed */
6873e197d65SBarry Smith       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
6883e197d65SBarry Smith       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
6893e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6903e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6913e197d65SBarry Smith       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
6923e197d65SBarry Smith       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6933e197d65SBarry Smith       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6943e197d65SBarry Smith     }
69551f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
69611755939SBarry Smith       cnt -= 2;
69751f519a2SBarry Smith       while (ilink->previous) {
69851f519a2SBarry Smith         ilink = ilink->previous;
69911755939SBarry Smith         /* compute the residual only over the part of the vector needed */
70011755939SBarry Smith         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
70111755939SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
70211755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
70311755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
70411755939SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
70511755939SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
70611755939SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
70751f519a2SBarry Smith       }
70811755939SBarry Smith     }
70965e19b50SBarry Smith   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
71051f519a2SBarry Smith   CHKMEMQ;
7110971522cSBarry Smith   PetscFunctionReturn(0);
7120971522cSBarry Smith }
7130971522cSBarry Smith 
714421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
715ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
716ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
717421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
718ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
719ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
720421e10b8SBarry Smith 
721421e10b8SBarry Smith #undef __FUNCT__
7228c03b21aSDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_FieldSplit"
723421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
724421e10b8SBarry Smith {
725421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
726421e10b8SBarry Smith   PetscErrorCode    ierr;
727421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
728421e10b8SBarry Smith 
729421e10b8SBarry Smith   PetscFunctionBegin;
730421e10b8SBarry Smith   CHKMEMQ;
731421e10b8SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
732421e10b8SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
733421e10b8SBarry Smith 
734421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
735421e10b8SBarry Smith     if (jac->defaultsplit) {
736421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
737421e10b8SBarry Smith       while (ilink) {
738421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
739421e10b8SBarry Smith 	ilink = ilink->next;
740421e10b8SBarry Smith       }
741421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
742421e10b8SBarry Smith     } else {
743421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
744421e10b8SBarry Smith       while (ilink) {
745421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
746421e10b8SBarry Smith 	ilink = ilink->next;
747421e10b8SBarry Smith       }
748421e10b8SBarry Smith     }
749421e10b8SBarry Smith   } else {
750421e10b8SBarry Smith     if (!jac->w1) {
751421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
752421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
753421e10b8SBarry Smith     }
754421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
755421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
756421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
757421e10b8SBarry Smith       while (ilink->next) {
758421e10b8SBarry Smith         ilink = ilink->next;
7599989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
760421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
761421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
762421e10b8SBarry Smith       }
763421e10b8SBarry Smith       while (ilink->previous) {
764421e10b8SBarry Smith         ilink = ilink->previous;
7659989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
766421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
767421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
768421e10b8SBarry Smith       }
769421e10b8SBarry Smith     } else {
770421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
771421e10b8SBarry Smith 	ilink = ilink->next;
772421e10b8SBarry Smith       }
773421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
774421e10b8SBarry Smith       while (ilink->previous) {
775421e10b8SBarry Smith 	ilink = ilink->previous;
7769989ab13SBarry Smith 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
777421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
778421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
779421e10b8SBarry Smith       }
780421e10b8SBarry Smith     }
781421e10b8SBarry Smith   }
782421e10b8SBarry Smith   CHKMEMQ;
783421e10b8SBarry Smith   PetscFunctionReturn(0);
784421e10b8SBarry Smith }
785421e10b8SBarry Smith 
7860971522cSBarry Smith #undef __FUNCT__
787574deadeSBarry Smith #define __FUNCT__ "PCReset_FieldSplit"
788574deadeSBarry Smith static PetscErrorCode PCReset_FieldSplit(PC pc)
7890971522cSBarry Smith {
7900971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
7910971522cSBarry Smith   PetscErrorCode    ierr;
7925a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
7930971522cSBarry Smith 
7940971522cSBarry Smith   PetscFunctionBegin;
7955a9f2f41SSatish Balay   while (ilink) {
796574deadeSBarry Smith     ierr = KSPReset(ilink->ksp);CHKERRQ(ierr);
797fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->x);CHKERRQ(ierr);
798fcfd50ebSBarry Smith     ierr = VecDestroy(&ilink->y);CHKERRQ(ierr);
799fcfd50ebSBarry Smith     ierr = VecScatterDestroy(&ilink->sctx);CHKERRQ(ierr);
800fcfd50ebSBarry Smith     ierr = ISDestroy(&ilink->is);CHKERRQ(ierr);
8015a9f2f41SSatish Balay     next = ilink->next;
8025a9f2f41SSatish Balay     ilink = next;
8030971522cSBarry Smith   }
80405b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
805574deadeSBarry Smith   if (jac->mat && jac->mat != jac->pmat) {
806574deadeSBarry Smith     ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);
807574deadeSBarry Smith   } else if (jac->mat) {
808574deadeSBarry Smith     jac->mat = PETSC_NULL;
809574deadeSBarry Smith   }
81097bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
81168dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
8126bf464f9SBarry Smith   ierr = VecDestroy(&jac->w1);CHKERRQ(ierr);
8136bf464f9SBarry Smith   ierr = VecDestroy(&jac->w2);CHKERRQ(ierr);
8146bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur);CHKERRQ(ierr);
8156bf464f9SBarry Smith   ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
8166bf464f9SBarry Smith   ierr = KSPDestroy(&jac->kspschur);CHKERRQ(ierr);
8176bf464f9SBarry Smith   ierr = MatDestroy(&jac->B);CHKERRQ(ierr);
8186bf464f9SBarry Smith   ierr = MatDestroy(&jac->C);CHKERRQ(ierr);
81963ec74ffSBarry Smith   jac->reset = PETSC_TRUE;
820574deadeSBarry Smith   PetscFunctionReturn(0);
821574deadeSBarry Smith }
822574deadeSBarry Smith 
823574deadeSBarry Smith #undef __FUNCT__
824574deadeSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
825574deadeSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
826574deadeSBarry Smith {
827574deadeSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
828574deadeSBarry Smith   PetscErrorCode    ierr;
829574deadeSBarry Smith   PC_FieldSplitLink ilink = jac->head,next;
830574deadeSBarry Smith 
831574deadeSBarry Smith   PetscFunctionBegin;
832574deadeSBarry Smith   ierr = PCReset_FieldSplit(pc);CHKERRQ(ierr);
833574deadeSBarry Smith   while (ilink) {
8346bf464f9SBarry Smith     ierr = KSPDestroy(&ilink->ksp);CHKERRQ(ierr);
835574deadeSBarry Smith     next = ilink->next;
836574deadeSBarry Smith     ierr = PetscFree(ilink->splitname);CHKERRQ(ierr);
837574deadeSBarry Smith     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
838574deadeSBarry Smith     ierr = PetscFree(ilink);CHKERRQ(ierr);
839574deadeSBarry Smith     ilink = next;
840574deadeSBarry Smith   }
841574deadeSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
842c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
843d88c663fSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","",PETSC_NULL);CHKERRQ(ierr);
844d88c663fSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","",PETSC_NULL);CHKERRQ(ierr);
845d88c663fSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","",PETSC_NULL);CHKERRQ(ierr);
846d88c663fSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","",PETSC_NULL);CHKERRQ(ierr);
847d88c663fSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","",PETSC_NULL);CHKERRQ(ierr);
8480971522cSBarry Smith   PetscFunctionReturn(0);
8490971522cSBarry Smith }
8500971522cSBarry Smith 
8510971522cSBarry Smith #undef __FUNCT__
8520971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
8530971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
8540971522cSBarry Smith {
8551b9fc7fcSBarry Smith   PetscErrorCode  ierr;
8566c924f48SJed Brown   PetscInt        bs;
857bc59fbc5SBarry Smith   PetscBool       flg,stokes = PETSC_FALSE;
8589dcbbd2bSBarry Smith   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
8593b224e63SBarry Smith   PCCompositeType ctype;
8601b9fc7fcSBarry Smith 
8610971522cSBarry Smith   PetscFunctionBegin;
8621b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
863acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr);
86451f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
86551f519a2SBarry Smith   if (flg) {
86651f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
86751f519a2SBarry Smith   }
868704ba839SBarry Smith 
869435f959eSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_fieldsplit_detect_saddle_point",&stokes,PETSC_NULL);CHKERRQ(ierr);
870c0adfefeSBarry Smith   if (stokes) {
871c0adfefeSBarry Smith     ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
872c0adfefeSBarry Smith     jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_SELF;
873c0adfefeSBarry Smith   }
874c0adfefeSBarry Smith 
8753b224e63SBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
8763b224e63SBarry Smith   if (flg) {
8773b224e63SBarry Smith     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
8783b224e63SBarry Smith   }
879d32f9abdSBarry Smith 
880c30613efSMatthew Knepley   /* Only setup fields once */
881c30613efSMatthew Knepley   if ((jac->bs > 0) && (jac->nsplits == 0)) {
882d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
883d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
8846c924f48SJed Brown     ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr);
8856c924f48SJed Brown     if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);}
886d32f9abdSBarry Smith   }
887c5d2311dSJed Brown   if (jac->type == PC_COMPOSITE_SCHUR) {
888c5d2311dSJed Brown     ierr = PetscOptionsEnum("-pc_fieldsplit_schur_factorization_type","Factorization to use","None",PCFieldSplitSchurFactorizationTypes,(PetscEnum)jac->schurfactorization,(PetscEnum*)&jac->schurfactorization,PETSC_NULL);CHKERRQ(ierr);
889084e4875SJed 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);
890c5d2311dSJed Brown   }
8911b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
8920971522cSBarry Smith   PetscFunctionReturn(0);
8930971522cSBarry Smith }
8940971522cSBarry Smith 
8950971522cSBarry Smith /*------------------------------------------------------------------------------------*/
8960971522cSBarry Smith 
8970971522cSBarry Smith EXTERN_C_BEGIN
8980971522cSBarry Smith #undef __FUNCT__
8990971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
9007087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields)
9010971522cSBarry Smith {
90297bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
9030971522cSBarry Smith   PetscErrorCode    ierr;
9045a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
90569a612a9SBarry Smith   char              prefix[128];
90651f519a2SBarry Smith   PetscInt          i;
9070971522cSBarry Smith 
9080971522cSBarry Smith   PetscFunctionBegin;
9096c924f48SJed Brown   if (jac->splitdefined) {
9106c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
9116c924f48SJed Brown     PetscFunctionReturn(0);
9126c924f48SJed Brown   }
91351f519a2SBarry Smith   for (i=0; i<n; i++) {
914e32f2f54SBarry 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);
915e32f2f54SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
91651f519a2SBarry Smith   }
917704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
918a04f6461SBarry Smith   if (splitname) {
919db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
920a04f6461SBarry Smith   } else {
921a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
922a04f6461SBarry Smith     ierr = PetscSNPrintf(ilink->splitname,2,"%s",jac->nsplits);CHKERRQ(ierr);
923a04f6461SBarry Smith   }
924704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
9255a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
9265a9f2f41SSatish Balay   ilink->nfields = n;
9275a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
9287adad957SLisandro Dalcin   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
9291cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
9305a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
9319005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
93269a612a9SBarry Smith 
933a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
9345a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
9350971522cSBarry Smith 
9360971522cSBarry Smith   if (!next) {
9375a9f2f41SSatish Balay     jac->head       = ilink;
93851f519a2SBarry Smith     ilink->previous = PETSC_NULL;
9390971522cSBarry Smith   } else {
9400971522cSBarry Smith     while (next->next) {
9410971522cSBarry Smith       next = next->next;
9420971522cSBarry Smith     }
9435a9f2f41SSatish Balay     next->next      = ilink;
94451f519a2SBarry Smith     ilink->previous = next;
9450971522cSBarry Smith   }
9460971522cSBarry Smith   jac->nsplits++;
9470971522cSBarry Smith   PetscFunctionReturn(0);
9480971522cSBarry Smith }
9490971522cSBarry Smith EXTERN_C_END
9500971522cSBarry Smith 
951e69d4d44SBarry Smith EXTERN_C_BEGIN
952e69d4d44SBarry Smith #undef __FUNCT__
953e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
9547087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
955e69d4d44SBarry Smith {
956e69d4d44SBarry Smith   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
957e69d4d44SBarry Smith   PetscErrorCode ierr;
958e69d4d44SBarry Smith 
959e69d4d44SBarry Smith   PetscFunctionBegin;
960519d70e2SJed Brown   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
961e69d4d44SBarry Smith   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
962e69d4d44SBarry Smith   (*subksp)[1] = jac->kspschur;
96313e0d083SBarry Smith   if (n) *n = jac->nsplits;
964e69d4d44SBarry Smith   PetscFunctionReturn(0);
965e69d4d44SBarry Smith }
966e69d4d44SBarry Smith EXTERN_C_END
9670971522cSBarry Smith 
9680971522cSBarry Smith EXTERN_C_BEGIN
9690971522cSBarry Smith #undef __FUNCT__
97069a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
9717087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
9720971522cSBarry Smith {
9730971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
9740971522cSBarry Smith   PetscErrorCode    ierr;
9750971522cSBarry Smith   PetscInt          cnt = 0;
9765a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
9770971522cSBarry Smith 
9780971522cSBarry Smith   PetscFunctionBegin;
9795d480477SMatthew G Knepley   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
9805a9f2f41SSatish Balay   while (ilink) {
9815a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
9825a9f2f41SSatish Balay     ilink = ilink->next;
9830971522cSBarry Smith   }
9845d480477SMatthew 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);
98513e0d083SBarry Smith   if (n) *n = jac->nsplits;
9860971522cSBarry Smith   PetscFunctionReturn(0);
9870971522cSBarry Smith }
9880971522cSBarry Smith EXTERN_C_END
9890971522cSBarry Smith 
990704ba839SBarry Smith EXTERN_C_BEGIN
991704ba839SBarry Smith #undef __FUNCT__
992704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
9937087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is)
994704ba839SBarry Smith {
995704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
996704ba839SBarry Smith   PetscErrorCode    ierr;
997704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
998704ba839SBarry Smith   char              prefix[128];
999704ba839SBarry Smith 
1000704ba839SBarry Smith   PetscFunctionBegin;
10016c924f48SJed Brown   if (jac->splitdefined) {
10026c924f48SJed Brown     ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr);
10036c924f48SJed Brown     PetscFunctionReturn(0);
10046c924f48SJed Brown   }
100516913363SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
1006a04f6461SBarry Smith   if (splitname) {
1007db4c96c1SJed Brown     ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr);
1008a04f6461SBarry Smith   } else {
1009a04f6461SBarry Smith     ierr = PetscMalloc(3*sizeof(char),&ilink->splitname);CHKERRQ(ierr);
1010d88c663fSJed Brown     ierr = PetscSNPrintf(ilink->splitname,3,"%D",jac->nsplits);CHKERRQ(ierr);
1011a04f6461SBarry Smith   }
10121ab39975SBarry Smith   ilink->is      = is;
10131ab39975SBarry Smith   ierr           = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
1014704ba839SBarry Smith   ilink->next    = PETSC_NULL;
1015704ba839SBarry Smith   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
10161cee3971SBarry Smith   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
1017704ba839SBarry Smith   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
10189005cf84SBarry Smith   ierr           = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->ksp);CHKERRQ(ierr);
1019704ba839SBarry Smith 
1020a04f6461SBarry Smith   ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr);
1021704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
1022704ba839SBarry Smith 
1023704ba839SBarry Smith   if (!next) {
1024704ba839SBarry Smith     jac->head       = ilink;
1025704ba839SBarry Smith     ilink->previous = PETSC_NULL;
1026704ba839SBarry Smith   } else {
1027704ba839SBarry Smith     while (next->next) {
1028704ba839SBarry Smith       next = next->next;
1029704ba839SBarry Smith     }
1030704ba839SBarry Smith     next->next      = ilink;
1031704ba839SBarry Smith     ilink->previous = next;
1032704ba839SBarry Smith   }
1033704ba839SBarry Smith   jac->nsplits++;
1034704ba839SBarry Smith 
1035704ba839SBarry Smith   PetscFunctionReturn(0);
1036704ba839SBarry Smith }
1037704ba839SBarry Smith EXTERN_C_END
1038704ba839SBarry Smith 
10390971522cSBarry Smith #undef __FUNCT__
10400971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
10410971522cSBarry Smith /*@
10420971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
10430971522cSBarry Smith 
1044ad4df100SBarry Smith     Logically Collective on PC
10450971522cSBarry Smith 
10460971522cSBarry Smith     Input Parameters:
10470971522cSBarry Smith +   pc  - the preconditioner context
1048a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
10490971522cSBarry Smith .   n - the number of fields in this split
1050db4c96c1SJed Brown -   fields - the fields in this split
10510971522cSBarry Smith 
10520971522cSBarry Smith     Level: intermediate
10530971522cSBarry Smith 
1054d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
1055d32f9abdSBarry Smith 
1056d32f9abdSBarry Smith      The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block
1057d32f9abdSBarry 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
1058d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
1059d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
1060d32f9abdSBarry Smith 
1061db4c96c1SJed Brown      This function is called once per split (it creates a new split each time).  Solve options
1062db4c96c1SJed Brown      for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1063db4c96c1SJed Brown 
1064d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
10650971522cSBarry Smith 
10660971522cSBarry Smith @*/
10677087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields)
10680971522cSBarry Smith {
10694ac538c5SBarry Smith   PetscErrorCode ierr;
10700971522cSBarry Smith 
10710971522cSBarry Smith   PetscFunctionBegin;
10720700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1073db4c96c1SJed Brown   PetscValidCharPointer(splitname,2);
1074db4c96c1SJed Brown   if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname);
1075db4c96c1SJed Brown   PetscValidIntPointer(fields,3);
10764ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetFields_C",(PC,const char[],PetscInt,const PetscInt *),(pc,splitname,n,fields));CHKERRQ(ierr);
10770971522cSBarry Smith   PetscFunctionReturn(0);
10780971522cSBarry Smith }
10790971522cSBarry Smith 
10800971522cSBarry Smith #undef __FUNCT__
1081704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
1082704ba839SBarry Smith /*@
1083704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
1084704ba839SBarry Smith 
1085ad4df100SBarry Smith     Logically Collective on PC
1086704ba839SBarry Smith 
1087704ba839SBarry Smith     Input Parameters:
1088704ba839SBarry Smith +   pc  - the preconditioner context
1089a04f6461SBarry Smith .   splitname - name of this split, if PETSC_NULL the number of the split is used
1090db4c96c1SJed Brown -   is - the index set that defines the vector elements in this field
1091704ba839SBarry Smith 
1092d32f9abdSBarry Smith 
1093a6ffb8dbSJed Brown     Notes:
1094a6ffb8dbSJed Brown     Use PCFieldSplitSetFields(), for fields defined by strided types.
1095a6ffb8dbSJed Brown 
1096db4c96c1SJed Brown     This function is called once per split (it creates a new split each time).  Solve options
1097db4c96c1SJed Brown     for this split will be available under the prefix -fieldsplit_SPLITNAME_.
1098d32f9abdSBarry Smith 
1099704ba839SBarry Smith     Level: intermediate
1100704ba839SBarry Smith 
1101704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
1102704ba839SBarry Smith 
1103704ba839SBarry Smith @*/
11047087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetIS(PC pc,const char splitname[],IS is)
1105704ba839SBarry Smith {
11064ac538c5SBarry Smith   PetscErrorCode ierr;
1107704ba839SBarry Smith 
1108704ba839SBarry Smith   PetscFunctionBegin;
11090700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1110d88c663fSJed Brown   if (splitname) PetscValidCharPointer(splitname,2);
1111db4c96c1SJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,3);
11124ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetIS_C",(PC,const char[],IS),(pc,splitname,is));CHKERRQ(ierr);
1113704ba839SBarry Smith   PetscFunctionReturn(0);
1114704ba839SBarry Smith }
1115704ba839SBarry Smith 
1116704ba839SBarry Smith #undef __FUNCT__
111757a9adfeSMatthew G Knepley #define __FUNCT__ "PCFieldSplitGetIS"
111857a9adfeSMatthew G Knepley /*@
111957a9adfeSMatthew G Knepley     PCFieldSplitGetIS - Retrieves the elements for a field as an IS
112057a9adfeSMatthew G Knepley 
112157a9adfeSMatthew G Knepley     Logically Collective on PC
112257a9adfeSMatthew G Knepley 
112357a9adfeSMatthew G Knepley     Input Parameters:
112457a9adfeSMatthew G Knepley +   pc  - the preconditioner context
112557a9adfeSMatthew G Knepley -   splitname - name of this split
112657a9adfeSMatthew G Knepley 
112757a9adfeSMatthew G Knepley     Output Parameter:
112857a9adfeSMatthew G Knepley -   is - the index set that defines the vector elements in this field, or PETSC_NULL if the field is not found
112957a9adfeSMatthew G Knepley 
113057a9adfeSMatthew G Knepley     Level: intermediate
113157a9adfeSMatthew G Knepley 
113257a9adfeSMatthew G Knepley .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetIS()
113357a9adfeSMatthew G Knepley 
113457a9adfeSMatthew G Knepley @*/
113557a9adfeSMatthew G Knepley PetscErrorCode PCFieldSplitGetIS(PC pc,const char splitname[],IS *is)
113657a9adfeSMatthew G Knepley {
113757a9adfeSMatthew G Knepley   PetscErrorCode ierr;
113857a9adfeSMatthew G Knepley 
113957a9adfeSMatthew G Knepley   PetscFunctionBegin;
114057a9adfeSMatthew G Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
114157a9adfeSMatthew G Knepley   PetscValidCharPointer(splitname,2);
114257a9adfeSMatthew G Knepley   PetscValidPointer(is,3);
114357a9adfeSMatthew G Knepley   {
114457a9adfeSMatthew G Knepley     PC_FieldSplit    *jac   = (PC_FieldSplit *) pc->data;
114557a9adfeSMatthew G Knepley     PC_FieldSplitLink ilink = jac->head;
114657a9adfeSMatthew G Knepley     PetscBool         found;
114757a9adfeSMatthew G Knepley 
114857a9adfeSMatthew G Knepley     *is = PETSC_NULL;
114957a9adfeSMatthew G Knepley     while(ilink) {
115057a9adfeSMatthew G Knepley       ierr = PetscStrcmp(ilink->splitname, splitname, &found);CHKERRQ(ierr);
115157a9adfeSMatthew G Knepley       if (found) {
115257a9adfeSMatthew G Knepley         *is = ilink->is;
115357a9adfeSMatthew G Knepley         break;
115457a9adfeSMatthew G Knepley       }
115557a9adfeSMatthew G Knepley       ilink = ilink->next;
115657a9adfeSMatthew G Knepley     }
115757a9adfeSMatthew G Knepley   }
115857a9adfeSMatthew G Knepley   PetscFunctionReturn(0);
115957a9adfeSMatthew G Knepley }
116057a9adfeSMatthew G Knepley 
116157a9adfeSMatthew G Knepley #undef __FUNCT__
116251f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
116351f519a2SBarry Smith /*@
116451f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
116551f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
116651f519a2SBarry Smith 
1167ad4df100SBarry Smith     Logically Collective on PC
116851f519a2SBarry Smith 
116951f519a2SBarry Smith     Input Parameters:
117051f519a2SBarry Smith +   pc  - the preconditioner context
117151f519a2SBarry Smith -   bs - the block size
117251f519a2SBarry Smith 
117351f519a2SBarry Smith     Level: intermediate
117451f519a2SBarry Smith 
117551f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
117651f519a2SBarry Smith 
117751f519a2SBarry Smith @*/
11787087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
117951f519a2SBarry Smith {
11804ac538c5SBarry Smith   PetscErrorCode ierr;
118151f519a2SBarry Smith 
118251f519a2SBarry Smith   PetscFunctionBegin;
11830700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1184c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,bs,2);
11854ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetBlockSize_C",(PC,PetscInt),(pc,bs));CHKERRQ(ierr);
118651f519a2SBarry Smith   PetscFunctionReturn(0);
118751f519a2SBarry Smith }
118851f519a2SBarry Smith 
118951f519a2SBarry Smith #undef __FUNCT__
119069a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
11910971522cSBarry Smith /*@C
119269a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
11930971522cSBarry Smith 
119469a612a9SBarry Smith    Collective on KSP
11950971522cSBarry Smith 
11960971522cSBarry Smith    Input Parameter:
11970971522cSBarry Smith .  pc - the preconditioner context
11980971522cSBarry Smith 
11990971522cSBarry Smith    Output Parameters:
120013e0d083SBarry Smith +  n - the number of splits
120169a612a9SBarry Smith -  pc - the array of KSP contexts
12020971522cSBarry Smith 
12030971522cSBarry Smith    Note:
1204d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
1205d32f9abdSBarry Smith    (not the KSP just the array that contains them).
12060971522cSBarry Smith 
120769a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
12080971522cSBarry Smith 
12090971522cSBarry Smith    Level: advanced
12100971522cSBarry Smith 
12110971522cSBarry Smith .seealso: PCFIELDSPLIT
12120971522cSBarry Smith @*/
12137087cfbeSBarry Smith PetscErrorCode  PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
12140971522cSBarry Smith {
12154ac538c5SBarry Smith   PetscErrorCode ierr;
12160971522cSBarry Smith 
12170971522cSBarry Smith   PetscFunctionBegin;
12180700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
121913e0d083SBarry Smith   if (n) PetscValidIntPointer(n,2);
12204ac538c5SBarry Smith   ierr = PetscUseMethod(pc,"PCFieldSplitGetSubKSP_C",(PC,PetscInt*,KSP **),(pc,n,subksp));CHKERRQ(ierr);
12210971522cSBarry Smith   PetscFunctionReturn(0);
12220971522cSBarry Smith }
12230971522cSBarry Smith 
1224e69d4d44SBarry Smith #undef __FUNCT__
1225e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition"
1226e69d4d44SBarry Smith /*@
1227e69d4d44SBarry Smith     PCFieldSplitSchurPrecondition -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
1228a04f6461SBarry Smith       A11 matrix. Otherwise no preconditioner is used.
1229e69d4d44SBarry Smith 
1230e69d4d44SBarry Smith     Collective on PC
1231e69d4d44SBarry Smith 
1232e69d4d44SBarry Smith     Input Parameters:
1233e69d4d44SBarry Smith +   pc  - the preconditioner context
1234fd1303e9SJungho Lee .   ptype - which matrix to use for preconditioning the Schur complement, PC_FIELDSPLIT_SCHUR_PRE_DIAG (diag) is default
1235084e4875SJed Brown -   userpre - matrix to use for preconditioning, or PETSC_NULL
1236084e4875SJed Brown 
1237e69d4d44SBarry Smith     Options Database:
1238084e4875SJed Brown .     -pc_fieldsplit_schur_precondition <self,user,diag> default is diag
1239e69d4d44SBarry Smith 
1240fd1303e9SJungho Lee     Notes:
1241fd1303e9SJungho Lee $    If ptype is
1242fd1303e9SJungho Lee $        user then the preconditioner for the Schur complement is generated by the provided matrix (pre argument
1243fd1303e9SJungho Lee $             to this function).
1244fd1303e9SJungho Lee $        diag then the preconditioner for the Schur complement is generated by the block diagonal part of the original
1245fd1303e9SJungho Lee $             matrix associated with the Schur complement (i.e. A11)
1246fd1303e9SJungho Lee $        self the preconditioner for the Schur complement is generated from the Schur complement matrix itself:
1247fd1303e9SJungho Lee $             The only preconditioner that currently works directly with the Schur complement matrix object is the PCLSC
1248fd1303e9SJungho Lee $             preconditioner
1249fd1303e9SJungho Lee 
1250fd1303e9SJungho Lee      When solving a saddle point problem, where the A11 block is identically zero, using diag as the ptype only makes sense
1251fd1303e9SJungho Lee     with the additional option -fieldsplit_1_pc_type none. Usually for saddle point problems one would use a ptype of self and
1252fd1303e9SJungho Lee     -fieldsplit_1_pc_type lsc which uses the least squares commutator compute a preconditioner for the Schur complement.
1253fd1303e9SJungho Lee 
1254fd1303e9SJungho Lee     Developer Notes: This is a terrible name, gives no good indication of what the function does and should also have Set in
1255fd1303e9SJungho Lee      the name since it sets a proceedure to use.
1256fd1303e9SJungho Lee 
1257e69d4d44SBarry Smith     Level: intermediate
1258e69d4d44SBarry Smith 
1259fd1303e9SJungho Lee .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType, PCLSC
1260e69d4d44SBarry Smith 
1261e69d4d44SBarry Smith @*/
12627087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1263e69d4d44SBarry Smith {
12644ac538c5SBarry Smith   PetscErrorCode ierr;
1265e69d4d44SBarry Smith 
1266e69d4d44SBarry Smith   PetscFunctionBegin;
12670700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12684ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSchurPrecondition_C",(PC,PCFieldSplitSchurPreType,Mat),(pc,ptype,pre));CHKERRQ(ierr);
1269e69d4d44SBarry Smith   PetscFunctionReturn(0);
1270e69d4d44SBarry Smith }
1271e69d4d44SBarry Smith 
1272e69d4d44SBarry Smith EXTERN_C_BEGIN
1273e69d4d44SBarry Smith #undef __FUNCT__
1274e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
12757087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1276e69d4d44SBarry Smith {
1277e69d4d44SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1278084e4875SJed Brown   PetscErrorCode  ierr;
1279e69d4d44SBarry Smith 
1280e69d4d44SBarry Smith   PetscFunctionBegin;
1281084e4875SJed Brown   jac->schurpre = ptype;
1282084e4875SJed Brown   if (pre) {
12836bf464f9SBarry Smith     ierr = MatDestroy(&jac->schur_user);CHKERRQ(ierr);
1284084e4875SJed Brown     jac->schur_user = pre;
1285084e4875SJed Brown     ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1286084e4875SJed Brown   }
1287e69d4d44SBarry Smith   PetscFunctionReturn(0);
1288e69d4d44SBarry Smith }
1289e69d4d44SBarry Smith EXTERN_C_END
1290e69d4d44SBarry Smith 
129130ad9308SMatthew Knepley #undef __FUNCT__
129230ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
129330ad9308SMatthew Knepley /*@C
12948c03b21aSDmitry Karpeev    PCFieldSplitGetSchurBlocks - Gets all matrix blocks for the Schur complement
129530ad9308SMatthew Knepley 
129630ad9308SMatthew Knepley    Collective on KSP
129730ad9308SMatthew Knepley 
129830ad9308SMatthew Knepley    Input Parameter:
129930ad9308SMatthew Knepley .  pc - the preconditioner context
130030ad9308SMatthew Knepley 
130130ad9308SMatthew Knepley    Output Parameters:
1302a04f6461SBarry Smith +  A00 - the (0,0) block
1303a04f6461SBarry Smith .  A01 - the (0,1) block
1304a04f6461SBarry Smith .  A10 - the (1,0) block
1305a04f6461SBarry Smith -  A11 - the (1,1) block
130630ad9308SMatthew Knepley 
130730ad9308SMatthew Knepley    Level: advanced
130830ad9308SMatthew Knepley 
130930ad9308SMatthew Knepley .seealso: PCFIELDSPLIT
131030ad9308SMatthew Knepley @*/
1311a04f6461SBarry Smith PetscErrorCode  PCFieldSplitGetSchurBlocks(PC pc,Mat *A00,Mat *A01,Mat *A10, Mat *A11)
131230ad9308SMatthew Knepley {
131330ad9308SMatthew Knepley   PC_FieldSplit *jac = (PC_FieldSplit *) pc->data;
131430ad9308SMatthew Knepley 
131530ad9308SMatthew Knepley   PetscFunctionBegin;
13160700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1317c1235816SBarry Smith   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
1318a04f6461SBarry Smith   if (A00) *A00 = jac->pmat[0];
1319a04f6461SBarry Smith   if (A01) *A01 = jac->B;
1320a04f6461SBarry Smith   if (A10) *A10 = jac->C;
1321a04f6461SBarry Smith   if (A11) *A11 = jac->pmat[1];
132230ad9308SMatthew Knepley   PetscFunctionReturn(0);
132330ad9308SMatthew Knepley }
132430ad9308SMatthew Knepley 
132579416396SBarry Smith EXTERN_C_BEGIN
132679416396SBarry Smith #undef __FUNCT__
132779416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
13287087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
132979416396SBarry Smith {
133079416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1331e69d4d44SBarry Smith   PetscErrorCode ierr;
133279416396SBarry Smith 
133379416396SBarry Smith   PetscFunctionBegin;
133479416396SBarry Smith   jac->type = type;
13353b224e63SBarry Smith   if (type == PC_COMPOSITE_SCHUR) {
13363b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit_Schur;
13373b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit_Schur;
1338e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1339e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1340e69d4d44SBarry Smith 
13413b224e63SBarry Smith   } else {
13423b224e63SBarry Smith     pc->ops->apply = PCApply_FieldSplit;
13433b224e63SBarry Smith     pc->ops->view  = PCView_FieldSplit;
1344e69d4d44SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
13459e7d6b0aSBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr);
13463b224e63SBarry Smith   }
134779416396SBarry Smith   PetscFunctionReturn(0);
134879416396SBarry Smith }
134979416396SBarry Smith EXTERN_C_END
135079416396SBarry Smith 
135151f519a2SBarry Smith EXTERN_C_BEGIN
135251f519a2SBarry Smith #undef __FUNCT__
135351f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
13547087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
135551f519a2SBarry Smith {
135651f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
135751f519a2SBarry Smith 
135851f519a2SBarry Smith   PetscFunctionBegin;
135965e19b50SBarry Smith   if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
136065e19b50SBarry 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);
136151f519a2SBarry Smith   jac->bs = bs;
136251f519a2SBarry Smith   PetscFunctionReturn(0);
136351f519a2SBarry Smith }
136451f519a2SBarry Smith EXTERN_C_END
136551f519a2SBarry Smith 
136679416396SBarry Smith #undef __FUNCT__
136779416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
1368bc08b0f1SBarry Smith /*@
136979416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
137079416396SBarry Smith 
137179416396SBarry Smith    Collective on PC
137279416396SBarry Smith 
137379416396SBarry Smith    Input Parameter:
137479416396SBarry Smith .  pc - the preconditioner context
137581540f2fSBarry Smith .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
137679416396SBarry Smith 
137779416396SBarry Smith    Options Database Key:
1378a4efd8eaSMatthew Knepley .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
137979416396SBarry Smith 
138079416396SBarry Smith    Level: Developer
138179416396SBarry Smith 
138279416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
138379416396SBarry Smith 
138479416396SBarry Smith .seealso: PCCompositeSetType()
138579416396SBarry Smith 
138679416396SBarry Smith @*/
13877087cfbeSBarry Smith PetscErrorCode  PCFieldSplitSetType(PC pc,PCCompositeType type)
138879416396SBarry Smith {
13894ac538c5SBarry Smith   PetscErrorCode ierr;
139079416396SBarry Smith 
139179416396SBarry Smith   PetscFunctionBegin;
13920700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13934ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFieldSplitSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr);
139479416396SBarry Smith   PetscFunctionReturn(0);
139579416396SBarry Smith }
139679416396SBarry Smith 
13970971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
13980971522cSBarry Smith /*MC
1399a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
1400a04f6461SBarry Smith                   fields or groups of fields. See the users manual section "Solving Block Matrices" for more details.
14010971522cSBarry Smith 
1402edf189efSBarry Smith      To set options on the solvers for each block append -fieldsplit_ to all the PC
1403edf189efSBarry Smith         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
14040971522cSBarry Smith 
1405a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
140669a612a9SBarry Smith          and set the options directly on the resulting KSP object
14070971522cSBarry Smith 
14080971522cSBarry Smith    Level: intermediate
14090971522cSBarry Smith 
141079416396SBarry Smith    Options Database Keys:
141181540f2fSBarry Smith +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
141281540f2fSBarry Smith .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
141381540f2fSBarry Smith                               been supplied explicitly by -pc_fieldsplit_%d_fields
141481540f2fSBarry Smith .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
141581540f2fSBarry Smith .   -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative>
1416e69d4d44SBarry Smith .   -pc_fieldsplit_schur_precondition <true,false> default is true
1417435f959eSBarry Smith .   -pc_fieldsplit_detect_saddle_point - automatically finds rows with zero or negative diagonal and uses Schur complement with no preconditioner as the solver
141843890ad2SJed Brown .   -fieldsplit_NAME_ksp_* - control inner linear solver, NAME is a sequential integer if unspecified, otherwise use name provided in PCFieldSplitSetIS() or the name associated with the field in the DM passed to PCSetDM() (or a higher level function)
141943890ad2SJed Brown -   -fieldsplit_NAME_pc_* - control inner preconditioner, NAME has same semantics as above
142079416396SBarry Smith 
142143890ad2SJed Brown    Notes:
142243890ad2SJed Brown     Use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1423d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
1424d32f9abdSBarry Smith 
1425d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
1426d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1427d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
1428d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
1429d32f9abdSBarry Smith 
1430b7bde157SMatthew G Knepley $     For the Schur complement preconditioner if J = ( A00 A01 )
1431b7bde157SMatthew G Knepley $                                                    ( A10 A11 )
1432b7bde157SMatthew G Knepley $     the preconditioner using full factorization is
1433b7bde157SMatthew G Knepley $              ( I   -A10 ksp(A00) ) ( inv(A00)     0  ) (     I          0  )
1434b7bde157SMatthew G Knepley $              ( 0         I       ) (   0      ksp(S) ) ( -A10 ksp(A00)  I  )
1435a04f6461SBarry Smith      where the action of inv(A00) is applied using the KSP solver with prefix -fieldsplit_0_. The action of
1436b7bde157SMatthew G Knepley      ksp(S) is computed using the KSP solver with prefix -fieldsplit_splitname_ (where splitname was given
1437b7bde157SMatthew G Knepley      in providing the SECOND split or 1 if not give). For PCFieldSplitGetKSP() when field number is 0,
1438b7bde157SMatthew G Knepley      it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
1439a04f6461SBarry Smith      A11 is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
1440b7bde157SMatthew G Knepley      option. You can use the preconditioner PCLSC to precondition the Schur complement with -fieldsplit_1_pc_type lsc. The
1441b7bde157SMatthew G Knepley      factorization type is set using -pc_fieldsplit_schur_factorization_type <diag, lower, upper, full>. The full is shown above,
1442b7bde157SMatthew G Knepley      but diag gives
1443b7bde157SMatthew G Knepley $              ( inv(A00)     0   )
1444b7bde157SMatthew G Knepley $              (   0      -ksp(S) )
1445b7bde157SMatthew G Knepley      so that the preconditioner is positive definite. The lower factorization is the inverse of
1446b7bde157SMatthew G Knepley $              (  A00   0 )
1447b7bde157SMatthew G Knepley $              (  A10   S )
1448b7bde157SMatthew G Knepley      where the inverses of A00 and S are applied using KSPs. The upper factorization is the inverse of
1449b7bde157SMatthew G Knepley $              ( A00 A01 )
1450b7bde157SMatthew G Knepley $              (  0   S  )
1451b7bde157SMatthew G Knepley      where again the inverses of A00 and S are applied using KSPs.
1452e69d4d44SBarry Smith 
1453edf189efSBarry Smith      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
1454edf189efSBarry Smith      is used automatically for a second block.
1455edf189efSBarry Smith 
1456ff218e97SBarry Smith      The fieldsplit preconditioner cannot currently be used with the BAIJ or SBAIJ data formats if the blocksize is larger than 1.
1457ff218e97SBarry Smith      Generally it should be used with the AIJ format.
1458ff218e97SBarry Smith 
1459ff218e97SBarry Smith      The forms of these preconditioners are closely related if not identical to forms derived as "Distributive Iterations", see,
1460ff218e97SBarry Smith      for example, page 294 in "Principles of Computational Fluid Dynamics" by Pieter Wesseling. Note that one can also use PCFIELDSPLIT
1461ff218e97SBarry Smith      inside a smoother resulting in "Distributive Smoothers".
14620716a85fSBarry Smith 
1463a541d17aSBarry Smith    Concepts: physics based preconditioners, block preconditioners
14640971522cSBarry Smith 
14657e8cb189SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCLSC,
1466e69d4d44SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
14670971522cSBarry Smith M*/
14680971522cSBarry Smith 
14690971522cSBarry Smith EXTERN_C_BEGIN
14700971522cSBarry Smith #undef __FUNCT__
14710971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
14727087cfbeSBarry Smith PetscErrorCode  PCCreate_FieldSplit(PC pc)
14730971522cSBarry Smith {
14740971522cSBarry Smith   PetscErrorCode ierr;
14750971522cSBarry Smith   PC_FieldSplit  *jac;
14760971522cSBarry Smith 
14770971522cSBarry Smith   PetscFunctionBegin;
147838f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
14790971522cSBarry Smith   jac->bs        = -1;
14800971522cSBarry Smith   jac->nsplits   = 0;
14813e197d65SBarry Smith   jac->type      = PC_COMPOSITE_MULTIPLICATIVE;
1482e6cab6aaSJed Brown   jac->schurpre  = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
1483c5d2311dSJed Brown   jac->schurfactorization = PC_FIELDSPLIT_SCHUR_FACTORIZATION_FULL;
148451f519a2SBarry Smith 
14850971522cSBarry Smith   pc->data     = (void*)jac;
14860971522cSBarry Smith 
14870971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
1488421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
14890971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
1490574deadeSBarry Smith   pc->ops->reset             = PCReset_FieldSplit;
14910971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
14920971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
14930971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
14940971522cSBarry Smith   pc->ops->applyrichardson   = 0;
14950971522cSBarry Smith 
149669a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
149769a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
14980971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
14990971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1500704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
1501704ba839SBarry Smith                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
150279416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
150379416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
150451f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
150551f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
15060971522cSBarry Smith   PetscFunctionReturn(0);
15070971522cSBarry Smith }
15080971522cSBarry Smith EXTERN_C_END
15090971522cSBarry Smith 
15100971522cSBarry Smith 
1511a541d17aSBarry Smith 
1512