1dba47a55SKris Buschelman #define PETSCKSP_DLL 2dba47a55SKris Buschelman 36356e834SBarry Smith #include "private/pcimpl.h" /*I "petscpc.h" I*/ 40971522cSBarry Smith 5084e4875SJed Brown const char *PCFieldSplitSchurPreTypes[] = {"self","diag","user","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0}; 6084e4875SJed Brown 70971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink; 80971522cSBarry Smith struct _PC_FieldSplitLink { 969a612a9SBarry Smith KSP ksp; 100971522cSBarry Smith Vec x,y; 110971522cSBarry Smith PetscInt nfields; 120971522cSBarry Smith PetscInt *fields; 131b9fc7fcSBarry Smith VecScatter sctx; 144aa3045dSJed Brown IS is; 1551f519a2SBarry Smith PC_FieldSplitLink next,previous; 160971522cSBarry Smith }; 170971522cSBarry Smith 180971522cSBarry Smith typedef struct { 1968dd23aaSBarry Smith PCCompositeType type; 20a4efd8eaSMatthew Knepley PetscTruth defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */ 21519d70e2SJed Brown PetscTruth realdiagonal; /* Flag to use the diagonal blocks of mat preconditioned by pmat, instead of just pmat */ 2230ad9308SMatthew Knepley PetscInt bs; /* Block size for IS and Mat structures */ 2330ad9308SMatthew Knepley PetscInt nsplits; /* Number of field divisions defined */ 2479416396SBarry Smith Vec *x,*y,w1,w2; 25519d70e2SJed Brown Mat *mat; /* The diagonal block for each split */ 26519d70e2SJed Brown Mat *pmat; /* The preconditioning diagonal block for each split */ 2730ad9308SMatthew Knepley Mat *Afield; /* The rows of the matrix associated with each split */ 28704ba839SBarry Smith PetscTruth issetup; 2930ad9308SMatthew Knepley /* Only used when Schur complement preconditioning is used */ 3030ad9308SMatthew Knepley Mat B; /* The (0,1) block */ 3130ad9308SMatthew Knepley Mat C; /* The (1,0) block */ 3230ad9308SMatthew Knepley Mat schur; /* The Schur complement S = D - C A^{-1} B */ 33084e4875SJed Brown Mat schur_user; /* User-provided preconditioning matrix for the Schur complement */ 34084e4875SJed Brown PCFieldSplitSchurPreType schurpre; /* Determines which preconditioning matrix is used for the Schur complement */ 3530ad9308SMatthew Knepley KSP kspschur; /* The solver for S */ 3697bbdb24SBarry Smith PC_FieldSplitLink head; 370971522cSBarry Smith } PC_FieldSplit; 380971522cSBarry Smith 3916913363SBarry Smith /* 4016913363SBarry Smith Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of 4116913363SBarry Smith inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the 4216913363SBarry Smith PC you could change this. 4316913363SBarry Smith */ 44084e4875SJed Brown 45e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it. This way the 46084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */ 47084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac) 48084e4875SJed Brown { 49084e4875SJed Brown switch (jac->schurpre) { 50084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur; 51084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1]; 52084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */ 53084e4875SJed Brown default: 54084e4875SJed Brown return jac->schur_user ? jac->schur_user : jac->pmat[1]; 55084e4875SJed Brown } 56084e4875SJed Brown } 57084e4875SJed Brown 58084e4875SJed Brown 590971522cSBarry Smith #undef __FUNCT__ 600971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit" 610971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 620971522cSBarry Smith { 630971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 640971522cSBarry Smith PetscErrorCode ierr; 650971522cSBarry Smith PetscTruth iascii; 660971522cSBarry Smith PetscInt i,j; 675a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 680971522cSBarry Smith 690971522cSBarry Smith PetscFunctionBegin; 700971522cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 710971522cSBarry Smith if (iascii) { 7251f519a2SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr); 7369a612a9SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 740971522cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 750971522cSBarry Smith for (i=0; i<jac->nsplits; i++) { 761ab39975SBarry Smith if (ilink->fields) { 770971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 7879416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 795a9f2f41SSatish Balay for (j=0; j<ilink->nfields; j++) { 8079416396SBarry Smith if (j > 0) { 8179416396SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 8279416396SBarry Smith } 835a9f2f41SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 840971522cSBarry Smith } 850971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 8679416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 871ab39975SBarry Smith } else { 881ab39975SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 891ab39975SBarry Smith } 905a9f2f41SSatish Balay ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 915a9f2f41SSatish Balay ilink = ilink->next; 920971522cSBarry Smith } 930971522cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 940971522cSBarry Smith } else { 95e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 960971522cSBarry Smith } 970971522cSBarry Smith PetscFunctionReturn(0); 980971522cSBarry Smith } 990971522cSBarry Smith 1000971522cSBarry Smith #undef __FUNCT__ 1013b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur" 1023b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer) 1033b224e63SBarry Smith { 1043b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1053b224e63SBarry Smith PetscErrorCode ierr; 1063b224e63SBarry Smith PetscTruth iascii; 1073b224e63SBarry Smith PetscInt i,j; 1083b224e63SBarry Smith PC_FieldSplitLink ilink = jac->head; 1093b224e63SBarry Smith KSP ksp; 1103b224e63SBarry Smith 1113b224e63SBarry Smith PetscFunctionBegin; 1123b224e63SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1133b224e63SBarry Smith if (iascii) { 1143b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, blocksize = %D\n",jac->bs);CHKERRQ(ierr); 1153b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Split info:\n");CHKERRQ(ierr); 1163b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1173b224e63SBarry Smith for (i=0; i<jac->nsplits; i++) { 1183b224e63SBarry Smith if (ilink->fields) { 1193b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 1203b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 1213b224e63SBarry Smith for (j=0; j<ilink->nfields; j++) { 1223b224e63SBarry Smith if (j > 0) { 1233b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 1243b224e63SBarry Smith } 1253b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 1263b224e63SBarry Smith } 1273b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1283b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1293b224e63SBarry Smith } else { 1303b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1313b224e63SBarry Smith } 1323b224e63SBarry Smith ilink = ilink->next; 1333b224e63SBarry Smith } 1343b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A block \n");CHKERRQ(ierr); 1353b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 13612cae6f2SJed Brown if (jac->schur) { 13712cae6f2SJed Brown ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 1383b224e63SBarry Smith ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 13912cae6f2SJed Brown } else { 14012cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 14112cae6f2SJed Brown } 1423b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1433b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = D - C inv(A) B \n");CHKERRQ(ierr); 1443b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 14512cae6f2SJed Brown if (jac->kspschur) { 1463b224e63SBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 14712cae6f2SJed Brown } else { 14812cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 14912cae6f2SJed Brown } 1503b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1513b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1523b224e63SBarry Smith } else { 153e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1543b224e63SBarry Smith } 1553b224e63SBarry Smith PetscFunctionReturn(0); 1563b224e63SBarry Smith } 1573b224e63SBarry Smith 1583b224e63SBarry Smith #undef __FUNCT__ 15969a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 16069a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 1610971522cSBarry Smith { 1620971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1630971522cSBarry Smith PetscErrorCode ierr; 1645a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 165d32f9abdSBarry Smith PetscInt i = 0,*ifields,nfields; 166d32f9abdSBarry Smith PetscTruth flg = PETSC_FALSE,*fields,flg2; 167d32f9abdSBarry Smith char optionname[128]; 1680971522cSBarry Smith 1690971522cSBarry Smith PetscFunctionBegin; 170d32f9abdSBarry Smith if (!ilink) { 171d32f9abdSBarry Smith 172521d7252SBarry Smith if (jac->bs <= 0) { 173704ba839SBarry Smith if (pc->pmat) { 174521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 175704ba839SBarry Smith } else { 176704ba839SBarry Smith jac->bs = 1; 177704ba839SBarry Smith } 178521d7252SBarry Smith } 179d32f9abdSBarry Smith 180d32f9abdSBarry Smith ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 181d32f9abdSBarry Smith if (!flg) { 182d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 183d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 184d32f9abdSBarry Smith flg = PETSC_TRUE; /* switched off automatically if user sets fields manually here */ 185d32f9abdSBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 186d32f9abdSBarry Smith while (PETSC_TRUE) { 187d32f9abdSBarry Smith sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 188d32f9abdSBarry Smith nfields = jac->bs; 189d32f9abdSBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg2);CHKERRQ(ierr); 190d32f9abdSBarry Smith if (!flg2) break; 191e32f2f54SBarry Smith if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 192d32f9abdSBarry Smith flg = PETSC_FALSE; 193d32f9abdSBarry Smith ierr = PCFieldSplitSetFields(pc,nfields,ifields);CHKERRQ(ierr); 194d32f9abdSBarry Smith } 195d32f9abdSBarry Smith ierr = PetscFree(ifields);CHKERRQ(ierr); 196d32f9abdSBarry Smith } 197d32f9abdSBarry Smith 198d32f9abdSBarry Smith if (flg) { 199d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 20079416396SBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr); 20179416396SBarry Smith ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr); 2025a9f2f41SSatish Balay while (ilink) { 2035a9f2f41SSatish Balay for (i=0; i<ilink->nfields; i++) { 2045a9f2f41SSatish Balay fields[ilink->fields[i]] = PETSC_TRUE; 205521d7252SBarry Smith } 2065a9f2f41SSatish Balay ilink = ilink->next; 20779416396SBarry Smith } 20897bbdb24SBarry Smith jac->defaultsplit = PETSC_TRUE; 20979416396SBarry Smith for (i=0; i<jac->bs; i++) { 21079416396SBarry Smith if (!fields[i]) { 21179416396SBarry Smith ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr); 21279416396SBarry Smith } else { 21379416396SBarry Smith jac->defaultsplit = PETSC_FALSE; 21479416396SBarry Smith } 21579416396SBarry Smith } 21668dd23aaSBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 217521d7252SBarry Smith } 218edf189efSBarry Smith } else if (jac->nsplits == 1) { 219edf189efSBarry Smith if (ilink->is) { 220edf189efSBarry Smith IS is2; 221edf189efSBarry Smith PetscInt nmin,nmax; 222edf189efSBarry Smith 223edf189efSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 224edf189efSBarry Smith ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 225edf189efSBarry Smith ierr = PCFieldSplitSetIS(pc,is2);CHKERRQ(ierr); 226edf189efSBarry Smith ierr = ISDestroy(is2);CHKERRQ(ierr); 227*e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 228edf189efSBarry Smith } 229*e7e72b3dSBarry Smith if (jac->nsplits < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields"); 23069a612a9SBarry Smith PetscFunctionReturn(0); 23169a612a9SBarry Smith } 23269a612a9SBarry Smith 23369a612a9SBarry Smith 23469a612a9SBarry Smith #undef __FUNCT__ 23569a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 23669a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 23769a612a9SBarry Smith { 23869a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 23969a612a9SBarry Smith PetscErrorCode ierr; 2405a9f2f41SSatish Balay PC_FieldSplitLink ilink; 241cf502942SBarry Smith PetscInt i,nsplit,ccsize; 24269a612a9SBarry Smith MatStructure flag = pc->flag; 2436c8605c2SJed Brown PetscTruth sorted; 24469a612a9SBarry Smith 24569a612a9SBarry Smith PetscFunctionBegin; 24669a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 24797bbdb24SBarry Smith nsplit = jac->nsplits; 2485a9f2f41SSatish Balay ilink = jac->head; 24997bbdb24SBarry Smith 25097bbdb24SBarry Smith /* get the matrices for each split */ 251704ba839SBarry Smith if (!jac->issetup) { 2521b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 25397bbdb24SBarry Smith 254704ba839SBarry Smith jac->issetup = PETSC_TRUE; 255704ba839SBarry Smith 256704ba839SBarry Smith /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 25751f519a2SBarry Smith bs = jac->bs; 25897bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 259cf502942SBarry Smith ierr = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr); 2601b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 2611b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 2621b9fc7fcSBarry Smith if (jac->defaultsplit) { 263704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 264704ba839SBarry Smith } else if (!ilink->is) { 265ccb205f8SBarry Smith if (ilink->nfields > 1) { 2665a9f2f41SSatish Balay PetscInt *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields; 2675a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 2681b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 2691b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 2701b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 27197bbdb24SBarry Smith } 27297bbdb24SBarry Smith } 273704ba839SBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr); 274ccb205f8SBarry Smith ierr = PetscFree(ii);CHKERRQ(ierr); 275ccb205f8SBarry Smith } else { 276704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 277ccb205f8SBarry Smith } 2783e197d65SBarry Smith } 279704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 280e32f2f54SBarry Smith if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split"); 281704ba839SBarry Smith ilink = ilink->next; 2821b9fc7fcSBarry Smith } 2831b9fc7fcSBarry Smith } 2841b9fc7fcSBarry Smith 285704ba839SBarry Smith ilink = jac->head; 28697bbdb24SBarry Smith if (!jac->pmat) { 287cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 288cf502942SBarry Smith for (i=0; i<nsplit; i++) { 2894aa3045dSJed Brown ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 290704ba839SBarry Smith ilink = ilink->next; 291cf502942SBarry Smith } 29297bbdb24SBarry Smith } else { 293cf502942SBarry Smith for (i=0; i<nsplit; i++) { 2944aa3045dSJed Brown ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 295704ba839SBarry Smith ilink = ilink->next; 296cf502942SBarry Smith } 29797bbdb24SBarry Smith } 298519d70e2SJed Brown if (jac->realdiagonal) { 299519d70e2SJed Brown ilink = jac->head; 300519d70e2SJed Brown if (!jac->mat) { 301519d70e2SJed Brown ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr); 302519d70e2SJed Brown for (i=0; i<nsplit; i++) { 303519d70e2SJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 304519d70e2SJed Brown ilink = ilink->next; 305519d70e2SJed Brown } 306519d70e2SJed Brown } else { 307519d70e2SJed Brown for (i=0; i<nsplit; i++) { 308519d70e2SJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 309519d70e2SJed Brown ilink = ilink->next; 310519d70e2SJed Brown } 311519d70e2SJed Brown } 312519d70e2SJed Brown } else { 313519d70e2SJed Brown jac->mat = jac->pmat; 314519d70e2SJed Brown } 31597bbdb24SBarry Smith 3166c8605c2SJed Brown if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 31768dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 31868dd23aaSBarry Smith ilink = jac->head; 31968dd23aaSBarry Smith if (!jac->Afield) { 32068dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 32168dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 3224aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 32368dd23aaSBarry Smith ilink = ilink->next; 32468dd23aaSBarry Smith } 32568dd23aaSBarry Smith } else { 32668dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 3274aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 32868dd23aaSBarry Smith ilink = ilink->next; 32968dd23aaSBarry Smith } 33068dd23aaSBarry Smith } 33168dd23aaSBarry Smith } 33268dd23aaSBarry Smith 3333b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 3343b224e63SBarry Smith IS ccis; 3354aa3045dSJed Brown PetscInt rstart,rend; 336*e7e72b3dSBarry Smith if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 33768dd23aaSBarry Smith 338e6cab6aaSJed Brown /* When extracting off-diagonal submatrices, we take complements from this range */ 339e6cab6aaSJed Brown ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 340e6cab6aaSJed Brown 3413b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 3423b224e63SBarry Smith if (jac->schur) { 3433b224e63SBarry Smith ilink = jac->head; 3444aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 3454aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 3463b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 3473b224e63SBarry Smith ilink = ilink->next; 3484aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 3494aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 3503b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 351519d70e2SJed Brown ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr); 352084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 3533b224e63SBarry Smith 3543b224e63SBarry Smith } else { 3551cee3971SBarry Smith KSP ksp; 3563b224e63SBarry Smith 3573b224e63SBarry Smith /* extract the B and C matrices */ 3583b224e63SBarry Smith ilink = jac->head; 3594aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 3604aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 3613b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 3623b224e63SBarry Smith ilink = ilink->next; 3634aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 3644aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 3653b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 366084e4875SJed Brown /* Better would be to use 'mat[0]' (diagonal block of the real matrix) preconditioned by pmat[0] */ 367519d70e2SJed Brown ierr = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr); 3681cee3971SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 369e69d4d44SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr); 3703b224e63SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 3713b224e63SBarry Smith 3723b224e63SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr); 3731cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 374084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 375084e4875SJed Brown if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 376084e4875SJed Brown PC pc; 377cd5f4a64SJed Brown ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr); 378084e4875SJed Brown ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 379084e4875SJed Brown /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 380e69d4d44SBarry Smith } 3813b224e63SBarry Smith ierr = KSPSetOptionsPrefix(jac->kspschur,((PetscObject)pc)->prefix);CHKERRQ(ierr); 382edf189efSBarry Smith ierr = KSPAppendOptionsPrefix(jac->kspschur,"fieldsplit_1_");CHKERRQ(ierr); 3833b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 3843b224e63SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 3853b224e63SBarry Smith 3863b224e63SBarry Smith ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr); 3873b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr); 3883b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr); 3893b224e63SBarry Smith ilink = jac->head; 3903b224e63SBarry Smith ilink->x = jac->x[0]; ilink->y = jac->y[0]; 3913b224e63SBarry Smith ilink = ilink->next; 3923b224e63SBarry Smith ilink->x = jac->x[1]; ilink->y = jac->y[1]; 3933b224e63SBarry Smith } 3943b224e63SBarry Smith } else { 39597bbdb24SBarry Smith /* set up the individual PCs */ 39697bbdb24SBarry Smith i = 0; 3975a9f2f41SSatish Balay ilink = jac->head; 3985a9f2f41SSatish Balay while (ilink) { 399519d70e2SJed Brown ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr); 4003b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 4015a9f2f41SSatish Balay ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr); 4025a9f2f41SSatish Balay ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 40397bbdb24SBarry Smith i++; 4045a9f2f41SSatish Balay ilink = ilink->next; 4050971522cSBarry Smith } 40697bbdb24SBarry Smith 40797bbdb24SBarry Smith /* create work vectors for each split */ 4081b9fc7fcSBarry Smith if (!jac->x) { 40997bbdb24SBarry Smith ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 4105a9f2f41SSatish Balay ilink = jac->head; 41197bbdb24SBarry Smith for (i=0; i<nsplit; i++) { 412906ed7ccSBarry Smith Vec *vl,*vr; 413906ed7ccSBarry Smith 414906ed7ccSBarry Smith ierr = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr); 415906ed7ccSBarry Smith ilink->x = *vr; 416906ed7ccSBarry Smith ilink->y = *vl; 417906ed7ccSBarry Smith ierr = PetscFree(vr);CHKERRQ(ierr); 418906ed7ccSBarry Smith ierr = PetscFree(vl);CHKERRQ(ierr); 4195a9f2f41SSatish Balay jac->x[i] = ilink->x; 4205a9f2f41SSatish Balay jac->y[i] = ilink->y; 4215a9f2f41SSatish Balay ilink = ilink->next; 42297bbdb24SBarry Smith } 4233b224e63SBarry Smith } 4243b224e63SBarry Smith } 4253b224e63SBarry Smith 4263b224e63SBarry Smith 427d0f46423SBarry Smith if (!jac->head->sctx) { 4283b224e63SBarry Smith Vec xtmp; 4293b224e63SBarry Smith 43079416396SBarry Smith /* compute scatter contexts needed by multiplicative versions and non-default splits */ 4311b9fc7fcSBarry Smith 4325a9f2f41SSatish Balay ilink = jac->head; 4331b9fc7fcSBarry Smith ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 4341b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 435704ba839SBarry Smith ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 4365a9f2f41SSatish Balay ilink = ilink->next; 4371b9fc7fcSBarry Smith } 4381b9fc7fcSBarry Smith ierr = VecDestroy(xtmp);CHKERRQ(ierr); 4391b9fc7fcSBarry Smith } 4400971522cSBarry Smith PetscFunctionReturn(0); 4410971522cSBarry Smith } 4420971522cSBarry Smith 4435a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 444ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 445ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 4465a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 447ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 448ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 44979416396SBarry Smith 4500971522cSBarry Smith #undef __FUNCT__ 4513b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 4523b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 4533b224e63SBarry Smith { 4543b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 4553b224e63SBarry Smith PetscErrorCode ierr; 4563b224e63SBarry Smith KSP ksp; 4573b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 4583b224e63SBarry Smith 4593b224e63SBarry Smith PetscFunctionBegin; 4603b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 4613b224e63SBarry Smith 4623b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4633b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4643b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 4653b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 4663b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 4673b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4683b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4693b224e63SBarry Smith 4703b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 4713b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4723b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4733b224e63SBarry Smith 4743b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 4753b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 4763b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 4773b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4783b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4793b224e63SBarry Smith 4803b224e63SBarry Smith PetscFunctionReturn(0); 4813b224e63SBarry Smith } 4823b224e63SBarry Smith 4833b224e63SBarry Smith #undef __FUNCT__ 4840971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 4850971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 4860971522cSBarry Smith { 4870971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 4880971522cSBarry Smith PetscErrorCode ierr; 4895a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 490af93645bSJed Brown PetscInt cnt; 4910971522cSBarry Smith 4920971522cSBarry Smith PetscFunctionBegin; 49351f519a2SBarry Smith CHKMEMQ; 49451f519a2SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 49551f519a2SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 49651f519a2SBarry Smith 49779416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 4981b9fc7fcSBarry Smith if (jac->defaultsplit) { 4990971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 5005a9f2f41SSatish Balay while (ilink) { 5015a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 5025a9f2f41SSatish Balay ilink = ilink->next; 5030971522cSBarry Smith } 5040971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 5051b9fc7fcSBarry Smith } else { 506efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 5075a9f2f41SSatish Balay while (ilink) { 5085a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 5095a9f2f41SSatish Balay ilink = ilink->next; 5101b9fc7fcSBarry Smith } 5111b9fc7fcSBarry Smith } 51216913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 51379416396SBarry Smith if (!jac->w1) { 51479416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 51579416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 51679416396SBarry Smith } 517efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 5185a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 5193e197d65SBarry Smith cnt = 1; 5205a9f2f41SSatish Balay while (ilink->next) { 5215a9f2f41SSatish Balay ilink = ilink->next; 5223e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 5233e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 5243e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 5253e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5263e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5273e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 5283e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5293e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5303e197d65SBarry Smith } 53151f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 53211755939SBarry Smith cnt -= 2; 53351f519a2SBarry Smith while (ilink->previous) { 53451f519a2SBarry Smith ilink = ilink->previous; 53511755939SBarry Smith /* compute the residual only over the part of the vector needed */ 53611755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 53711755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 53811755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 53911755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 54011755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 54111755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 54211755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 54351f519a2SBarry Smith } 54411755939SBarry Smith } 545e32f2f54SBarry Smith } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 54651f519a2SBarry Smith CHKMEMQ; 5470971522cSBarry Smith PetscFunctionReturn(0); 5480971522cSBarry Smith } 5490971522cSBarry Smith 550421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 551ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 552ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 553421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 554ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 555ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 556421e10b8SBarry Smith 557421e10b8SBarry Smith #undef __FUNCT__ 558421e10b8SBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 559421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 560421e10b8SBarry Smith { 561421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 562421e10b8SBarry Smith PetscErrorCode ierr; 563421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 564421e10b8SBarry Smith 565421e10b8SBarry Smith PetscFunctionBegin; 566421e10b8SBarry Smith CHKMEMQ; 567421e10b8SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 568421e10b8SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 569421e10b8SBarry Smith 570421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 571421e10b8SBarry Smith if (jac->defaultsplit) { 572421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 573421e10b8SBarry Smith while (ilink) { 574421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 575421e10b8SBarry Smith ilink = ilink->next; 576421e10b8SBarry Smith } 577421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 578421e10b8SBarry Smith } else { 579421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 580421e10b8SBarry Smith while (ilink) { 581421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 582421e10b8SBarry Smith ilink = ilink->next; 583421e10b8SBarry Smith } 584421e10b8SBarry Smith } 585421e10b8SBarry Smith } else { 586421e10b8SBarry Smith if (!jac->w1) { 587421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 588421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 589421e10b8SBarry Smith } 590421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 591421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 592421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 593421e10b8SBarry Smith while (ilink->next) { 594421e10b8SBarry Smith ilink = ilink->next; 5959989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 596421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 597421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 598421e10b8SBarry Smith } 599421e10b8SBarry Smith while (ilink->previous) { 600421e10b8SBarry Smith ilink = ilink->previous; 6019989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 602421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 603421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 604421e10b8SBarry Smith } 605421e10b8SBarry Smith } else { 606421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 607421e10b8SBarry Smith ilink = ilink->next; 608421e10b8SBarry Smith } 609421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 610421e10b8SBarry Smith while (ilink->previous) { 611421e10b8SBarry Smith ilink = ilink->previous; 6129989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 613421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 614421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 615421e10b8SBarry Smith } 616421e10b8SBarry Smith } 617421e10b8SBarry Smith } 618421e10b8SBarry Smith CHKMEMQ; 619421e10b8SBarry Smith PetscFunctionReturn(0); 620421e10b8SBarry Smith } 621421e10b8SBarry Smith 6220971522cSBarry Smith #undef __FUNCT__ 6230971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 6240971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 6250971522cSBarry Smith { 6260971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6270971522cSBarry Smith PetscErrorCode ierr; 6285a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 6290971522cSBarry Smith 6300971522cSBarry Smith PetscFunctionBegin; 6315a9f2f41SSatish Balay while (ilink) { 6325a9f2f41SSatish Balay ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr); 6335a9f2f41SSatish Balay if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);} 6345a9f2f41SSatish Balay if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);} 6355a9f2f41SSatish Balay if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);} 636704ba839SBarry Smith if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);} 6375a9f2f41SSatish Balay next = ilink->next; 638704ba839SBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 639704ba839SBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 6405a9f2f41SSatish Balay ilink = next; 6410971522cSBarry Smith } 64205b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 643519d70e2SJed Brown if (jac->mat && jac->mat != jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);} 64497bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 64568dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 64679416396SBarry Smith if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);} 64779416396SBarry Smith if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);} 6483b224e63SBarry Smith if (jac->schur) {ierr = MatDestroy(jac->schur);CHKERRQ(ierr);} 649084e4875SJed Brown if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);} 650d0f46423SBarry Smith if (jac->kspschur) {ierr = KSPDestroy(jac->kspschur);CHKERRQ(ierr);} 6513b224e63SBarry Smith if (jac->B) {ierr = MatDestroy(jac->B);CHKERRQ(ierr);} 6523b224e63SBarry Smith if (jac->C) {ierr = MatDestroy(jac->C);CHKERRQ(ierr);} 6530971522cSBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 6540971522cSBarry Smith PetscFunctionReturn(0); 6550971522cSBarry Smith } 6560971522cSBarry Smith 6570971522cSBarry Smith #undef __FUNCT__ 6580971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 6590971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 6600971522cSBarry Smith { 6611b9fc7fcSBarry Smith PetscErrorCode ierr; 66251f519a2SBarry Smith PetscInt i = 0,nfields,*fields,bs; 663084e4875SJed Brown PetscTruth flg; 6641b9fc7fcSBarry Smith char optionname[128]; 6659dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6663b224e63SBarry Smith PCCompositeType ctype; 6671b9fc7fcSBarry Smith 6680971522cSBarry Smith PetscFunctionBegin; 6691b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 670519d70e2SJed Brown ierr = PetscOptionsTruth("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr); 67151f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 67251f519a2SBarry Smith if (flg) { 67351f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 67451f519a2SBarry Smith } 675704ba839SBarry Smith 6763b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 6773b224e63SBarry Smith if (flg) { 6783b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 6793b224e63SBarry Smith } 680d32f9abdSBarry Smith 681c30613efSMatthew Knepley /* Only setup fields once */ 682c30613efSMatthew Knepley if ((jac->bs > 0) && (jac->nsplits == 0)) { 683d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 684d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 68551f519a2SBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr); 6861b9fc7fcSBarry Smith while (PETSC_TRUE) { 68713f74950SBarry Smith sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 68851f519a2SBarry Smith nfields = jac->bs; 6891b9fc7fcSBarry Smith ierr = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr); 6901b9fc7fcSBarry Smith if (!flg) break; 691e32f2f54SBarry Smith if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 6921b9fc7fcSBarry Smith ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr); 6931b9fc7fcSBarry Smith } 69451f519a2SBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 695d32f9abdSBarry Smith } 696084e4875SJed 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); 6971b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 6980971522cSBarry Smith PetscFunctionReturn(0); 6990971522cSBarry Smith } 7000971522cSBarry Smith 7010971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 7020971522cSBarry Smith 7030971522cSBarry Smith EXTERN_C_BEGIN 7040971522cSBarry Smith #undef __FUNCT__ 7050971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 706dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields) 7070971522cSBarry Smith { 70897bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7090971522cSBarry Smith PetscErrorCode ierr; 7105a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 71169a612a9SBarry Smith char prefix[128]; 71251f519a2SBarry Smith PetscInt i; 7130971522cSBarry Smith 7140971522cSBarry Smith PetscFunctionBegin; 715*e7e72b3dSBarry Smith if (n <= 0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested"); 71651f519a2SBarry Smith for (i=0; i<n; i++) { 717e32f2f54SBarry 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); 718e32f2f54SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 71951f519a2SBarry Smith } 720704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 721704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 7225a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 7235a9f2f41SSatish Balay ilink->nfields = n; 7245a9f2f41SSatish Balay ilink->next = PETSC_NULL; 7257adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 7261cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 7275a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 72869a612a9SBarry Smith 7297adad957SLisandro Dalcin if (((PetscObject)pc)->prefix) { 7307adad957SLisandro Dalcin sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 73169a612a9SBarry Smith } else { 73213f74950SBarry Smith sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 73369a612a9SBarry Smith } 7345a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 7350971522cSBarry Smith 7360971522cSBarry Smith if (!next) { 7375a9f2f41SSatish Balay jac->head = ilink; 73851f519a2SBarry Smith ilink->previous = PETSC_NULL; 7390971522cSBarry Smith } else { 7400971522cSBarry Smith while (next->next) { 7410971522cSBarry Smith next = next->next; 7420971522cSBarry Smith } 7435a9f2f41SSatish Balay next->next = ilink; 74451f519a2SBarry Smith ilink->previous = next; 7450971522cSBarry Smith } 7460971522cSBarry Smith jac->nsplits++; 7470971522cSBarry Smith PetscFunctionReturn(0); 7480971522cSBarry Smith } 7490971522cSBarry Smith EXTERN_C_END 7500971522cSBarry Smith 751e69d4d44SBarry Smith EXTERN_C_BEGIN 752e69d4d44SBarry Smith #undef __FUNCT__ 753e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 754e69d4d44SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 755e69d4d44SBarry Smith { 756e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 757e69d4d44SBarry Smith PetscErrorCode ierr; 758e69d4d44SBarry Smith 759e69d4d44SBarry Smith PetscFunctionBegin; 760519d70e2SJed Brown ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 761e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 762e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 763084e4875SJed Brown *n = jac->nsplits; 764e69d4d44SBarry Smith PetscFunctionReturn(0); 765e69d4d44SBarry Smith } 766e69d4d44SBarry Smith EXTERN_C_END 7670971522cSBarry Smith 7680971522cSBarry Smith EXTERN_C_BEGIN 7690971522cSBarry Smith #undef __FUNCT__ 77069a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 771dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 7720971522cSBarry Smith { 7730971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7740971522cSBarry Smith PetscErrorCode ierr; 7750971522cSBarry Smith PetscInt cnt = 0; 7765a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 7770971522cSBarry Smith 7780971522cSBarry Smith PetscFunctionBegin; 77969a612a9SBarry Smith ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 7805a9f2f41SSatish Balay while (ilink) { 7815a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 7825a9f2f41SSatish Balay ilink = ilink->next; 7830971522cSBarry Smith } 784e32f2f54SBarry Smith if (cnt != jac->nsplits) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits); 7850971522cSBarry Smith *n = jac->nsplits; 7860971522cSBarry Smith PetscFunctionReturn(0); 7870971522cSBarry Smith } 7880971522cSBarry Smith EXTERN_C_END 7890971522cSBarry Smith 790704ba839SBarry Smith EXTERN_C_BEGIN 791704ba839SBarry Smith #undef __FUNCT__ 792704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 793704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,IS is) 794704ba839SBarry Smith { 795704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 796704ba839SBarry Smith PetscErrorCode ierr; 797704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 798704ba839SBarry Smith char prefix[128]; 799704ba839SBarry Smith 800704ba839SBarry Smith PetscFunctionBegin; 80116913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 8021ab39975SBarry Smith ilink->is = is; 8031ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 804704ba839SBarry Smith ilink->next = PETSC_NULL; 805704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 8061cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 807704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 808704ba839SBarry Smith 809704ba839SBarry Smith if (((PetscObject)pc)->prefix) { 810704ba839SBarry Smith sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 811704ba839SBarry Smith } else { 812704ba839SBarry Smith sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 813704ba839SBarry Smith } 814704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 815704ba839SBarry Smith 816704ba839SBarry Smith if (!next) { 817704ba839SBarry Smith jac->head = ilink; 818704ba839SBarry Smith ilink->previous = PETSC_NULL; 819704ba839SBarry Smith } else { 820704ba839SBarry Smith while (next->next) { 821704ba839SBarry Smith next = next->next; 822704ba839SBarry Smith } 823704ba839SBarry Smith next->next = ilink; 824704ba839SBarry Smith ilink->previous = next; 825704ba839SBarry Smith } 826704ba839SBarry Smith jac->nsplits++; 827704ba839SBarry Smith 828704ba839SBarry Smith PetscFunctionReturn(0); 829704ba839SBarry Smith } 830704ba839SBarry Smith EXTERN_C_END 831704ba839SBarry Smith 8320971522cSBarry Smith #undef __FUNCT__ 8330971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 8340971522cSBarry Smith /*@ 8350971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 8360971522cSBarry Smith 8370971522cSBarry Smith Collective on PC 8380971522cSBarry Smith 8390971522cSBarry Smith Input Parameters: 8400971522cSBarry Smith + pc - the preconditioner context 8410971522cSBarry Smith . n - the number of fields in this split 8420971522cSBarry Smith . fields - the fields in this split 8430971522cSBarry Smith 8440971522cSBarry Smith Level: intermediate 8450971522cSBarry Smith 846d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 847d32f9abdSBarry Smith 848d32f9abdSBarry Smith The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block 849d32f9abdSBarry 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 850d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 851d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 852d32f9abdSBarry Smith 853d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 8540971522cSBarry Smith 8550971522cSBarry Smith @*/ 856dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields) 8570971522cSBarry Smith { 8580971522cSBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *); 8590971522cSBarry Smith 8600971522cSBarry Smith PetscFunctionBegin; 8610700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 8620971522cSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr); 8630971522cSBarry Smith if (f) { 8640971522cSBarry Smith ierr = (*f)(pc,n,fields);CHKERRQ(ierr); 8650971522cSBarry Smith } 8660971522cSBarry Smith PetscFunctionReturn(0); 8670971522cSBarry Smith } 8680971522cSBarry Smith 8690971522cSBarry Smith #undef __FUNCT__ 870704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 871704ba839SBarry Smith /*@ 872704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 873704ba839SBarry Smith 874704ba839SBarry Smith Collective on PC 875704ba839SBarry Smith 876704ba839SBarry Smith Input Parameters: 877704ba839SBarry Smith + pc - the preconditioner context 878704ba839SBarry Smith . is - the index set that defines the vector elements in this field 879704ba839SBarry Smith 880d32f9abdSBarry Smith 881a6ffb8dbSJed Brown Notes: 882a6ffb8dbSJed Brown Use PCFieldSplitSetFields(), for fields defined by strided types. 883a6ffb8dbSJed Brown 884a6ffb8dbSJed Brown This function is called once per split (it creates a new split each time). 885d32f9abdSBarry Smith 886704ba839SBarry Smith Level: intermediate 887704ba839SBarry Smith 888704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 889704ba839SBarry Smith 890704ba839SBarry Smith @*/ 891704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,IS is) 892704ba839SBarry Smith { 893704ba839SBarry Smith PetscErrorCode ierr,(*f)(PC,IS); 894704ba839SBarry Smith 895704ba839SBarry Smith PetscFunctionBegin; 8960700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 8970700a824SBarry Smith PetscValidHeaderSpecific(is,IS_CLASSID,1); 898704ba839SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr); 899704ba839SBarry Smith if (f) { 900704ba839SBarry Smith ierr = (*f)(pc,is);CHKERRQ(ierr); 901704ba839SBarry Smith } 902704ba839SBarry Smith PetscFunctionReturn(0); 903704ba839SBarry Smith } 904704ba839SBarry Smith 905704ba839SBarry Smith #undef __FUNCT__ 90651f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 90751f519a2SBarry Smith /*@ 90851f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 90951f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 91051f519a2SBarry Smith 91151f519a2SBarry Smith Collective on PC 91251f519a2SBarry Smith 91351f519a2SBarry Smith Input Parameters: 91451f519a2SBarry Smith + pc - the preconditioner context 91551f519a2SBarry Smith - bs - the block size 91651f519a2SBarry Smith 91751f519a2SBarry Smith Level: intermediate 91851f519a2SBarry Smith 91951f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 92051f519a2SBarry Smith 92151f519a2SBarry Smith @*/ 92251f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 92351f519a2SBarry Smith { 92451f519a2SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt); 92551f519a2SBarry Smith 92651f519a2SBarry Smith PetscFunctionBegin; 9270700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 92851f519a2SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr); 92951f519a2SBarry Smith if (f) { 93051f519a2SBarry Smith ierr = (*f)(pc,bs);CHKERRQ(ierr); 93151f519a2SBarry Smith } 93251f519a2SBarry Smith PetscFunctionReturn(0); 93351f519a2SBarry Smith } 93451f519a2SBarry Smith 93551f519a2SBarry Smith #undef __FUNCT__ 93669a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 9370971522cSBarry Smith /*@C 93869a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 9390971522cSBarry Smith 94069a612a9SBarry Smith Collective on KSP 9410971522cSBarry Smith 9420971522cSBarry Smith Input Parameter: 9430971522cSBarry Smith . pc - the preconditioner context 9440971522cSBarry Smith 9450971522cSBarry Smith Output Parameters: 9460971522cSBarry Smith + n - the number of split 94769a612a9SBarry Smith - pc - the array of KSP contexts 9480971522cSBarry Smith 9490971522cSBarry Smith Note: 950d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 951d32f9abdSBarry Smith (not the KSP just the array that contains them). 9520971522cSBarry Smith 95369a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 9540971522cSBarry Smith 9550971522cSBarry Smith Level: advanced 9560971522cSBarry Smith 9570971522cSBarry Smith .seealso: PCFIELDSPLIT 9580971522cSBarry Smith @*/ 959dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 9600971522cSBarry Smith { 96169a612a9SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **); 9620971522cSBarry Smith 9630971522cSBarry Smith PetscFunctionBegin; 9640700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9650971522cSBarry Smith PetscValidIntPointer(n,2); 96669a612a9SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 9670971522cSBarry Smith if (f) { 96869a612a9SBarry Smith ierr = (*f)(pc,n,subksp);CHKERRQ(ierr); 969*e7e72b3dSBarry Smith } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC"); 9700971522cSBarry Smith PetscFunctionReturn(0); 9710971522cSBarry Smith } 9720971522cSBarry Smith 973e69d4d44SBarry Smith #undef __FUNCT__ 974e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 975e69d4d44SBarry Smith /*@ 976e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 977e69d4d44SBarry Smith D matrix. Otherwise no preconditioner is used. 978e69d4d44SBarry Smith 979e69d4d44SBarry Smith Collective on PC 980e69d4d44SBarry Smith 981e69d4d44SBarry Smith Input Parameters: 982e69d4d44SBarry Smith + pc - the preconditioner context 983084e4875SJed Brown . ptype - which matrix to use for preconditioning the Schur complement 984084e4875SJed Brown - userpre - matrix to use for preconditioning, or PETSC_NULL 985084e4875SJed Brown 986084e4875SJed Brown Notes: 987084e4875SJed Brown The default is to use the block on the diagonal of the preconditioning matrix. This is D, in the (1,1) position. 988084e4875SJed Brown There are currently no preconditioners that work directly with the Schur complement so setting 989084e4875SJed Brown PC_FIELDSPLIT_SCHUR_PRE_SELF is observationally equivalent to -fieldsplit_1_pc_type none. 990e69d4d44SBarry Smith 991e69d4d44SBarry Smith Options Database: 992084e4875SJed Brown . -pc_fieldsplit_schur_precondition <self,user,diag> default is diag 993e69d4d44SBarry Smith 994e69d4d44SBarry Smith Level: intermediate 995e69d4d44SBarry Smith 996084e4875SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType 997e69d4d44SBarry Smith 998e69d4d44SBarry Smith @*/ 999084e4875SJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1000e69d4d44SBarry Smith { 1001084e4875SJed Brown PetscErrorCode ierr,(*f)(PC,PCFieldSplitSchurPreType,Mat); 1002e69d4d44SBarry Smith 1003e69d4d44SBarry Smith PetscFunctionBegin; 10040700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1005e69d4d44SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",(void (**)(void))&f);CHKERRQ(ierr); 1006e69d4d44SBarry Smith if (f) { 1007084e4875SJed Brown ierr = (*f)(pc,ptype,pre);CHKERRQ(ierr); 1008e69d4d44SBarry Smith } 1009e69d4d44SBarry Smith PetscFunctionReturn(0); 1010e69d4d44SBarry Smith } 1011e69d4d44SBarry Smith 1012e69d4d44SBarry Smith EXTERN_C_BEGIN 1013e69d4d44SBarry Smith #undef __FUNCT__ 1014e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 1015084e4875SJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1016e69d4d44SBarry Smith { 1017e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1018084e4875SJed Brown PetscErrorCode ierr; 1019e69d4d44SBarry Smith 1020e69d4d44SBarry Smith PetscFunctionBegin; 1021084e4875SJed Brown jac->schurpre = ptype; 1022084e4875SJed Brown if (pre) { 1023084e4875SJed Brown if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);} 1024084e4875SJed Brown jac->schur_user = pre; 1025084e4875SJed Brown ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1026084e4875SJed Brown } 1027e69d4d44SBarry Smith PetscFunctionReturn(0); 1028e69d4d44SBarry Smith } 1029e69d4d44SBarry Smith EXTERN_C_END 1030e69d4d44SBarry Smith 103130ad9308SMatthew Knepley #undef __FUNCT__ 103230ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 103330ad9308SMatthew Knepley /*@C 103430ad9308SMatthew Knepley PCFieldSplitGetSchurBlocks - Gets the all matrix blocks for the Schur complement 103530ad9308SMatthew Knepley 103630ad9308SMatthew Knepley Collective on KSP 103730ad9308SMatthew Knepley 103830ad9308SMatthew Knepley Input Parameter: 103930ad9308SMatthew Knepley . pc - the preconditioner context 104030ad9308SMatthew Knepley 104130ad9308SMatthew Knepley Output Parameters: 104230ad9308SMatthew Knepley + A - the (0,0) block 104330ad9308SMatthew Knepley . B - the (0,1) block 104430ad9308SMatthew Knepley . C - the (1,0) block 104530ad9308SMatthew Knepley - D - the (1,1) block 104630ad9308SMatthew Knepley 104730ad9308SMatthew Knepley Level: advanced 104830ad9308SMatthew Knepley 104930ad9308SMatthew Knepley .seealso: PCFIELDSPLIT 105030ad9308SMatthew Knepley @*/ 105130ad9308SMatthew Knepley PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSchurBlocks(PC pc,Mat *A,Mat *B,Mat *C, Mat *D) 105230ad9308SMatthew Knepley { 105330ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 105430ad9308SMatthew Knepley 105530ad9308SMatthew Knepley PetscFunctionBegin; 10560700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1057*e7e72b3dSBarry Smith if (jac->type != PC_COMPOSITE_SCHUR) {SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");} 105830ad9308SMatthew Knepley if (A) *A = jac->pmat[0]; 105930ad9308SMatthew Knepley if (B) *B = jac->B; 106030ad9308SMatthew Knepley if (C) *C = jac->C; 106130ad9308SMatthew Knepley if (D) *D = jac->pmat[1]; 106230ad9308SMatthew Knepley PetscFunctionReturn(0); 106330ad9308SMatthew Knepley } 106430ad9308SMatthew Knepley 106579416396SBarry Smith EXTERN_C_BEGIN 106679416396SBarry Smith #undef __FUNCT__ 106779416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 1068dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 106979416396SBarry Smith { 107079416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1071e69d4d44SBarry Smith PetscErrorCode ierr; 107279416396SBarry Smith 107379416396SBarry Smith PetscFunctionBegin; 107479416396SBarry Smith jac->type = type; 10753b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 10763b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 10773b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 1078e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1079e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1080e69d4d44SBarry Smith 10813b224e63SBarry Smith } else { 10823b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 10833b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 1084e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 10859e7d6b0aSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 10863b224e63SBarry Smith } 108779416396SBarry Smith PetscFunctionReturn(0); 108879416396SBarry Smith } 108979416396SBarry Smith EXTERN_C_END 109079416396SBarry Smith 109151f519a2SBarry Smith EXTERN_C_BEGIN 109251f519a2SBarry Smith #undef __FUNCT__ 109351f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 109451f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 109551f519a2SBarry Smith { 109651f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 109751f519a2SBarry Smith 109851f519a2SBarry Smith PetscFunctionBegin; 1099e32f2f54SBarry Smith if (bs < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 1100e32f2f54SBarry Smith if (jac->bs > 0 && jac->bs != bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Cannot change fieldsplit blocksize from %D to %D after it has been set",jac->bs,bs); 110151f519a2SBarry Smith jac->bs = bs; 110251f519a2SBarry Smith PetscFunctionReturn(0); 110351f519a2SBarry Smith } 110451f519a2SBarry Smith EXTERN_C_END 110551f519a2SBarry Smith 110679416396SBarry Smith #undef __FUNCT__ 110779416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1108bc08b0f1SBarry Smith /*@ 110979416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 111079416396SBarry Smith 111179416396SBarry Smith Collective on PC 111279416396SBarry Smith 111379416396SBarry Smith Input Parameter: 111479416396SBarry Smith . pc - the preconditioner context 111581540f2fSBarry Smith . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 111679416396SBarry Smith 111779416396SBarry Smith Options Database Key: 1118a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 111979416396SBarry Smith 112079416396SBarry Smith Level: Developer 112179416396SBarry Smith 112279416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 112379416396SBarry Smith 112479416396SBarry Smith .seealso: PCCompositeSetType() 112579416396SBarry Smith 112679416396SBarry Smith @*/ 1127dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type) 112879416396SBarry Smith { 112979416396SBarry Smith PetscErrorCode ierr,(*f)(PC,PCCompositeType); 113079416396SBarry Smith 113179416396SBarry Smith PetscFunctionBegin; 11320700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 113379416396SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 113479416396SBarry Smith if (f) { 113579416396SBarry Smith ierr = (*f)(pc,type);CHKERRQ(ierr); 113679416396SBarry Smith } 113779416396SBarry Smith PetscFunctionReturn(0); 113879416396SBarry Smith } 113979416396SBarry Smith 11400971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 11410971522cSBarry Smith /*MC 1142a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 11430971522cSBarry Smith fields or groups of fields 11440971522cSBarry Smith 11450971522cSBarry Smith 1146edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1147edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 11480971522cSBarry Smith 1149a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 115069a612a9SBarry Smith and set the options directly on the resulting KSP object 11510971522cSBarry Smith 11520971522cSBarry Smith Level: intermediate 11530971522cSBarry Smith 115479416396SBarry Smith Options Database Keys: 115581540f2fSBarry Smith + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 115681540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 115781540f2fSBarry Smith been supplied explicitly by -pc_fieldsplit_%d_fields 115881540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 115981540f2fSBarry Smith . -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative> 1160e69d4d44SBarry Smith . -pc_fieldsplit_schur_precondition <true,false> default is true 116179416396SBarry Smith 1162edf189efSBarry Smith - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 11633b224e63SBarry Smith for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 11643b224e63SBarry Smith 11653b224e63SBarry Smith 1166d32f9abdSBarry Smith Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1167d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1168d32f9abdSBarry Smith 1169d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1170d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1171d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1172d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1173d32f9abdSBarry Smith 1174d32f9abdSBarry Smith Currently for the multiplicative version, the updated residual needed for the next field 1175d32f9abdSBarry Smith solve is computed via a matrix vector product over the entire array. An optimization would be 1176d32f9abdSBarry Smith to update the residual only for the part of the right hand side associated with the next field 1177d32f9abdSBarry Smith solve. (This would involve more MatGetSubMatrix() calls or some other mechanism to compute the 1178d32f9abdSBarry Smith part of the matrix needed to just update part of the residual). 1179d32f9abdSBarry Smith 1180e69d4d44SBarry Smith For the Schur complement preconditioner if J = ( A B ) 1181e69d4d44SBarry Smith ( C D ) 1182e69d4d44SBarry Smith the preconditioner is 1183e69d4d44SBarry Smith (I -B inv(A)) ( inv(A) 0 ) (I 0 ) 1184e69d4d44SBarry Smith (0 I ) ( 0 inv(S) ) (-C inv(A) I ) 1185edf189efSBarry Smith where the action of inv(A) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1186e69d4d44SBarry Smith inv(S) is computed using the KSP solver with prefix -schur_. For PCFieldSplitGetKSP() when field number is 1187edf189efSBarry Smith 0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1188e69d4d44SBarry Smith D is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1189e69d4d44SBarry Smith option. 1190e69d4d44SBarry Smith 1191edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1192edf189efSBarry Smith is used automatically for a second block. 1193edf189efSBarry Smith 1194a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 11950971522cSBarry Smith 1196a541d17aSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners 1197e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 11980971522cSBarry Smith M*/ 11990971522cSBarry Smith 12000971522cSBarry Smith EXTERN_C_BEGIN 12010971522cSBarry Smith #undef __FUNCT__ 12020971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 1203dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc) 12040971522cSBarry Smith { 12050971522cSBarry Smith PetscErrorCode ierr; 12060971522cSBarry Smith PC_FieldSplit *jac; 12070971522cSBarry Smith 12080971522cSBarry Smith PetscFunctionBegin; 120938f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 12100971522cSBarry Smith jac->bs = -1; 12110971522cSBarry Smith jac->nsplits = 0; 12123e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1213e6cab6aaSJed Brown jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 121451f519a2SBarry Smith 12150971522cSBarry Smith pc->data = (void*)jac; 12160971522cSBarry Smith 12170971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1218421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 12190971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 12200971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 12210971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 12220971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 12230971522cSBarry Smith pc->ops->applyrichardson = 0; 12240971522cSBarry Smith 122569a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 122669a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 12270971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 12280971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1229704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1230704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 123179416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 123279416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 123351f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 123451f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 12350971522cSBarry Smith PetscFunctionReturn(0); 12360971522cSBarry Smith } 12370971522cSBarry Smith EXTERN_C_END 12380971522cSBarry Smith 12390971522cSBarry Smith 1240a541d17aSBarry Smith /*MC 1241a541d17aSBarry Smith Block_Preconditioners - PETSc provides support for a variety of "block preconditioners", this provides an 1242a541d17aSBarry Smith overview of these methods. 1243a541d17aSBarry Smith 1244a541d17aSBarry Smith Consider the solution to ( A B ) (x_1) = (b_1) 1245a541d17aSBarry Smith ( C D ) (x_2) (b_2) 1246a541d17aSBarry Smith 1247a541d17aSBarry Smith Important special cases, the Stokes equation: C = B' and D = 0 (A B) (x_1) = (b_1) 1248a541d17aSBarry Smith B' 0) (x_2) (b_2) 1249a541d17aSBarry Smith 1250a541d17aSBarry Smith One of the goals of the PCFieldSplit preconditioner in PETSc is to provide a variety of preconditioners 1251a541d17aSBarry Smith for this block system. 1252a541d17aSBarry Smith 1253a541d17aSBarry Smith Consider an additional matrix (Ap Bp) 1254a541d17aSBarry Smith (Cp Dp) where some or all of the entries may be the same as 1255a541d17aSBarry Smith in the original matrix (for example Ap == A). 1256a541d17aSBarry Smith 1257a541d17aSBarry Smith In the following, A^ denotes the approximate application of the inverse of A, possibly using Ap in the 1258a541d17aSBarry Smith approximation. In PETSc this simply means one has called KSPSetOperators(ksp,A,Ap,...) or KSPSetOperators(ksp,Ap,Ap,...) 1259a541d17aSBarry Smith 1260a541d17aSBarry Smith Block Jacobi: x_1 = A^ b_1 1261a541d17aSBarry Smith x_2 = D^ b_2 1262a541d17aSBarry Smith 1263a541d17aSBarry Smith Lower block Gauss-Seidel: x_1 = A^ b_1 1264a541d17aSBarry Smith x_2 = D^ (b_2 - C x_1) variant x_2 = D^ (b_2 - Cp x_1) 1265a541d17aSBarry Smith 1266a541d17aSBarry Smith Symmetric Gauss-Seidel: x_1 = x_1 + A^(b_1 - A x_1 - B x_2) variant x_1 = x_1 + A^(b_1 - Ap x_1 - Bp x_2) 1267a541d17aSBarry Smith Interestingly this form is not actually a symmetric matrix, the symmetric version is 1268a541d17aSBarry Smith x_1 = A^(b_1 - B x_2) variant x_1 = A^(b_1 - Bp x_2) 1269a541d17aSBarry Smith 12700bc0a719SBarry Smith Level: intermediate 12710bc0a719SBarry Smith 1272a541d17aSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCFIELDSPLIT 1273a541d17aSBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 1274a541d17aSBarry Smith M*/ 1275