1dba47a55SKris Buschelman #define PETSCKSP_DLL 2dba47a55SKris Buschelman 30971522cSBarry Smith /* 40971522cSBarry Smith 50971522cSBarry Smith */ 66356e834SBarry Smith #include "private/pcimpl.h" /*I "petscpc.h" I*/ 70971522cSBarry Smith 8084e4875SJed Brown const char *PCFieldSplitSchurPreTypes[] = {"self","diag","user","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0}; 9084e4875SJed Brown 100971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink; 110971522cSBarry Smith struct _PC_FieldSplitLink { 1269a612a9SBarry Smith KSP ksp; 130971522cSBarry Smith Vec x,y; 140971522cSBarry Smith PetscInt nfields; 150971522cSBarry Smith PetscInt *fields; 161b9fc7fcSBarry Smith VecScatter sctx; 174aa3045dSJed Brown IS is; 1851f519a2SBarry Smith PC_FieldSplitLink next,previous; 190971522cSBarry Smith }; 200971522cSBarry Smith 210971522cSBarry Smith typedef struct { 2268dd23aaSBarry Smith PCCompositeType type; 23a4efd8eaSMatthew Knepley PetscTruth defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */ 2430ad9308SMatthew Knepley PetscInt bs; /* Block size for IS and Mat structures */ 2530ad9308SMatthew Knepley PetscInt nsplits; /* Number of field divisions defined */ 2679416396SBarry Smith Vec *x,*y,w1,w2; 2730ad9308SMatthew Knepley Mat *pmat; /* The diagonal block for each split */ 2830ad9308SMatthew Knepley Mat *Afield; /* The rows of the matrix associated with each split */ 29704ba839SBarry Smith PetscTruth issetup; 3030ad9308SMatthew Knepley /* Only used when Schur complement preconditioning is used */ 3130ad9308SMatthew Knepley Mat B; /* The (0,1) block */ 3230ad9308SMatthew Knepley Mat C; /* The (1,0) block */ 3330ad9308SMatthew Knepley Mat schur; /* The Schur complement S = D - C A^{-1} B */ 34084e4875SJed Brown Mat schur_user; /* User-provided preconditioning matrix for the Schur complement */ 35084e4875SJed Brown PCFieldSplitSchurPreType schurpre; /* Determines which preconditioning matrix is used for the Schur complement */ 3630ad9308SMatthew Knepley KSP kspschur; /* The solver for S */ 3797bbdb24SBarry Smith PC_FieldSplitLink head; 380971522cSBarry Smith } PC_FieldSplit; 390971522cSBarry Smith 4016913363SBarry Smith /* 4116913363SBarry Smith Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of 4216913363SBarry Smith inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the 4316913363SBarry Smith PC you could change this. 4416913363SBarry Smith */ 45084e4875SJed Brown 46*e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it. This way the 47084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */ 48084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac) 49084e4875SJed Brown { 50084e4875SJed Brown switch (jac->schurpre) { 51084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur; 52084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1]; 53084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */ 54084e4875SJed Brown default: 55084e4875SJed Brown return jac->schur_user ? jac->schur_user : jac->pmat[1]; 56084e4875SJed Brown } 57084e4875SJed Brown } 58084e4875SJed Brown 59084e4875SJed Brown 600971522cSBarry Smith #undef __FUNCT__ 610971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit" 620971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 630971522cSBarry Smith { 640971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 650971522cSBarry Smith PetscErrorCode ierr; 660971522cSBarry Smith PetscTruth iascii; 670971522cSBarry Smith PetscInt i,j; 685a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 690971522cSBarry Smith 700971522cSBarry Smith PetscFunctionBegin; 710971522cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 720971522cSBarry Smith if (iascii) { 7351f519a2SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr); 7469a612a9SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 750971522cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 760971522cSBarry Smith for (i=0; i<jac->nsplits; i++) { 771ab39975SBarry Smith if (ilink->fields) { 780971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 7979416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 805a9f2f41SSatish Balay for (j=0; j<ilink->nfields; j++) { 8179416396SBarry Smith if (j > 0) { 8279416396SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 8379416396SBarry Smith } 845a9f2f41SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 850971522cSBarry Smith } 860971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 8779416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 881ab39975SBarry Smith } else { 891ab39975SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 901ab39975SBarry Smith } 915a9f2f41SSatish Balay ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 925a9f2f41SSatish Balay ilink = ilink->next; 930971522cSBarry Smith } 940971522cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 950971522cSBarry Smith } else { 960971522cSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 970971522cSBarry Smith } 980971522cSBarry Smith PetscFunctionReturn(0); 990971522cSBarry Smith } 1000971522cSBarry Smith 1010971522cSBarry Smith #undef __FUNCT__ 1023b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur" 1033b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer) 1043b224e63SBarry Smith { 1053b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1063b224e63SBarry Smith PetscErrorCode ierr; 1073b224e63SBarry Smith PetscTruth iascii; 1083b224e63SBarry Smith PetscInt i,j; 1093b224e63SBarry Smith PC_FieldSplitLink ilink = jac->head; 1103b224e63SBarry Smith KSP ksp; 1113b224e63SBarry Smith 1123b224e63SBarry Smith PetscFunctionBegin; 1133b224e63SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1143b224e63SBarry Smith if (iascii) { 1153b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, blocksize = %D\n",jac->bs);CHKERRQ(ierr); 1163b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Split info:\n");CHKERRQ(ierr); 1173b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1183b224e63SBarry Smith for (i=0; i<jac->nsplits; i++) { 1193b224e63SBarry Smith if (ilink->fields) { 1203b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 1213b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 1223b224e63SBarry Smith for (j=0; j<ilink->nfields; j++) { 1233b224e63SBarry Smith if (j > 0) { 1243b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 1253b224e63SBarry Smith } 1263b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 1273b224e63SBarry Smith } 1283b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1293b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1303b224e63SBarry Smith } else { 1313b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1323b224e63SBarry Smith } 1333b224e63SBarry Smith ilink = ilink->next; 1343b224e63SBarry Smith } 1353b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A block \n");CHKERRQ(ierr); 1363b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 1373b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1383b224e63SBarry Smith ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 1393b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1403b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = D - C inv(A) B \n");CHKERRQ(ierr); 1413b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1423b224e63SBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 1433b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1443b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1453b224e63SBarry Smith } else { 1463b224e63SBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1473b224e63SBarry Smith } 1483b224e63SBarry Smith PetscFunctionReturn(0); 1493b224e63SBarry Smith } 1503b224e63SBarry Smith 1513b224e63SBarry Smith #undef __FUNCT__ 15269a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 15369a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 1540971522cSBarry Smith { 1550971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1560971522cSBarry Smith PetscErrorCode ierr; 1575a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 158d32f9abdSBarry Smith PetscInt i = 0,*ifields,nfields; 159d32f9abdSBarry Smith PetscTruth flg = PETSC_FALSE,*fields,flg2; 160d32f9abdSBarry Smith char optionname[128]; 1610971522cSBarry Smith 1620971522cSBarry Smith PetscFunctionBegin; 163d32f9abdSBarry Smith if (!ilink) { 164d32f9abdSBarry Smith 165521d7252SBarry Smith if (jac->bs <= 0) { 166704ba839SBarry Smith if (pc->pmat) { 167521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 168704ba839SBarry Smith } else { 169704ba839SBarry Smith jac->bs = 1; 170704ba839SBarry Smith } 171521d7252SBarry Smith } 172d32f9abdSBarry Smith 173d32f9abdSBarry Smith ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 174d32f9abdSBarry Smith if (!flg) { 175d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 176d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 177d32f9abdSBarry Smith flg = PETSC_TRUE; /* switched off automatically if user sets fields manually here */ 178d32f9abdSBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 179d32f9abdSBarry Smith while (PETSC_TRUE) { 180d32f9abdSBarry Smith sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 181d32f9abdSBarry Smith nfields = jac->bs; 182d32f9abdSBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg2);CHKERRQ(ierr); 183d32f9abdSBarry Smith if (!flg2) break; 184d32f9abdSBarry Smith if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 185d32f9abdSBarry Smith flg = PETSC_FALSE; 186d32f9abdSBarry Smith ierr = PCFieldSplitSetFields(pc,nfields,ifields);CHKERRQ(ierr); 187d32f9abdSBarry Smith } 188d32f9abdSBarry Smith ierr = PetscFree(ifields);CHKERRQ(ierr); 189d32f9abdSBarry Smith } 190d32f9abdSBarry Smith 191d32f9abdSBarry Smith if (flg) { 192d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 19379416396SBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr); 19479416396SBarry Smith ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr); 1955a9f2f41SSatish Balay while (ilink) { 1965a9f2f41SSatish Balay for (i=0; i<ilink->nfields; i++) { 1975a9f2f41SSatish Balay fields[ilink->fields[i]] = PETSC_TRUE; 198521d7252SBarry Smith } 1995a9f2f41SSatish Balay ilink = ilink->next; 20079416396SBarry Smith } 20197bbdb24SBarry Smith jac->defaultsplit = PETSC_TRUE; 20279416396SBarry Smith for (i=0; i<jac->bs; i++) { 20379416396SBarry Smith if (!fields[i]) { 20479416396SBarry Smith ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr); 20579416396SBarry Smith } else { 20679416396SBarry Smith jac->defaultsplit = PETSC_FALSE; 20779416396SBarry Smith } 20879416396SBarry Smith } 20968dd23aaSBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 210521d7252SBarry Smith } 211edf189efSBarry Smith } else if (jac->nsplits == 1) { 212edf189efSBarry Smith if (ilink->is) { 213edf189efSBarry Smith IS is2; 214edf189efSBarry Smith PetscInt nmin,nmax; 215edf189efSBarry Smith 216edf189efSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 217edf189efSBarry Smith ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 218edf189efSBarry Smith ierr = PCFieldSplitSetIS(pc,is2);CHKERRQ(ierr); 219edf189efSBarry Smith ierr = ISDestroy(is2);CHKERRQ(ierr); 220edf189efSBarry Smith } else { 221edf189efSBarry Smith SETERRQ(PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 222edf189efSBarry Smith } 223d32f9abdSBarry Smith } 22469a612a9SBarry Smith PetscFunctionReturn(0); 22569a612a9SBarry Smith } 22669a612a9SBarry Smith 22769a612a9SBarry Smith 22869a612a9SBarry Smith #undef __FUNCT__ 22969a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 23069a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 23169a612a9SBarry Smith { 23269a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 23369a612a9SBarry Smith PetscErrorCode ierr; 2345a9f2f41SSatish Balay PC_FieldSplitLink ilink; 235cf502942SBarry Smith PetscInt i,nsplit,ccsize; 23669a612a9SBarry Smith MatStructure flag = pc->flag; 23768dd23aaSBarry Smith PetscTruth sorted,getsub; 23869a612a9SBarry Smith 23969a612a9SBarry Smith PetscFunctionBegin; 24069a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 24197bbdb24SBarry Smith nsplit = jac->nsplits; 2425a9f2f41SSatish Balay ilink = jac->head; 24397bbdb24SBarry Smith 24497bbdb24SBarry Smith /* get the matrices for each split */ 245704ba839SBarry Smith if (!jac->issetup) { 2461b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 24797bbdb24SBarry Smith 248704ba839SBarry Smith jac->issetup = PETSC_TRUE; 249704ba839SBarry Smith 250704ba839SBarry Smith /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 25151f519a2SBarry Smith bs = jac->bs; 25297bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 253cf502942SBarry Smith ierr = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr); 2541b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 2551b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 2561b9fc7fcSBarry Smith if (jac->defaultsplit) { 257704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 258704ba839SBarry Smith } else if (!ilink->is) { 259ccb205f8SBarry Smith if (ilink->nfields > 1) { 2605a9f2f41SSatish Balay PetscInt *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields; 2615a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 2621b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 2631b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 2641b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 26597bbdb24SBarry Smith } 26697bbdb24SBarry Smith } 267704ba839SBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr); 268ccb205f8SBarry Smith ierr = PetscFree(ii);CHKERRQ(ierr); 269ccb205f8SBarry Smith } else { 270704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 271ccb205f8SBarry Smith } 2723e197d65SBarry Smith } 273704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 2741b9fc7fcSBarry Smith if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split"); 275704ba839SBarry Smith ilink = ilink->next; 2761b9fc7fcSBarry Smith } 2771b9fc7fcSBarry Smith } 2781b9fc7fcSBarry Smith 279704ba839SBarry Smith ilink = jac->head; 28097bbdb24SBarry Smith if (!jac->pmat) { 281cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 282cf502942SBarry Smith for (i=0; i<nsplit; i++) { 2834aa3045dSJed Brown ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 284704ba839SBarry Smith ilink = ilink->next; 285cf502942SBarry Smith } 28697bbdb24SBarry Smith } else { 287cf502942SBarry Smith for (i=0; i<nsplit; i++) { 2884aa3045dSJed Brown ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 289704ba839SBarry Smith ilink = ilink->next; 290cf502942SBarry Smith } 29197bbdb24SBarry Smith } 29297bbdb24SBarry Smith 293*e6cab6aaSJed Brown if (jac->type != PC_COMPOSITE_SCHUR) { 29468dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 29568dd23aaSBarry Smith ierr = MatHasOperation(pc->mat,MATOP_GET_SUBMATRIX,&getsub);CHKERRQ(ierr); 2963b224e63SBarry Smith if (getsub && jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 29768dd23aaSBarry Smith ilink = jac->head; 29868dd23aaSBarry Smith if (!jac->Afield) { 29968dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 30068dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 3014aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 30268dd23aaSBarry Smith ilink = ilink->next; 30368dd23aaSBarry Smith } 30468dd23aaSBarry Smith } else { 30568dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 3064aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 30768dd23aaSBarry Smith ilink = ilink->next; 30868dd23aaSBarry Smith } 30968dd23aaSBarry Smith } 31068dd23aaSBarry Smith } 311*e6cab6aaSJed Brown } 31268dd23aaSBarry Smith 3133b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 3143b224e63SBarry Smith IS ccis; 3154aa3045dSJed Brown PetscInt rstart,rend; 3163b224e63SBarry Smith if (nsplit != 2) SETERRQ(PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 31768dd23aaSBarry Smith 318*e6cab6aaSJed Brown /* When extracting off-diagonal submatrices, we take complements from this range */ 319*e6cab6aaSJed Brown ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 320*e6cab6aaSJed Brown 3213b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 3223b224e63SBarry Smith if (jac->schur) { 3233b224e63SBarry Smith ilink = jac->head; 3244aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 3254aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 3263b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 3273b224e63SBarry Smith ilink = ilink->next; 3284aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 3294aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 3303b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 331084e4875SJed Brown ierr = MatSchurComplementUpdate(jac->schur,jac->pmat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr); 332084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 3333b224e63SBarry Smith 3343b224e63SBarry Smith } else { 3351cee3971SBarry Smith KSP ksp; 3363b224e63SBarry Smith 3373b224e63SBarry Smith /* extract the B and C matrices */ 3383b224e63SBarry Smith ilink = jac->head; 3394aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 3404aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 3413b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 3423b224e63SBarry Smith ilink = ilink->next; 3434aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 3444aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 3453b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 346084e4875SJed Brown /* Better would be to use 'mat[0]' (diagonal block of the real matrix) preconditioned by pmat[0] */ 347084e4875SJed Brown ierr = MatCreateSchurComplement(jac->pmat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],&jac->schur);CHKERRQ(ierr); 3481cee3971SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 349e69d4d44SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr); 3503b224e63SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 3513b224e63SBarry Smith 3523b224e63SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr); 3531cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 354084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 355084e4875SJed Brown if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 356084e4875SJed Brown PC pc; 357084e4875SJed Brown ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 358084e4875SJed Brown ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 359084e4875SJed Brown /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 360e69d4d44SBarry Smith } 3613b224e63SBarry Smith ierr = KSPSetOptionsPrefix(jac->kspschur,((PetscObject)pc)->prefix);CHKERRQ(ierr); 362edf189efSBarry Smith ierr = KSPAppendOptionsPrefix(jac->kspschur,"fieldsplit_1_");CHKERRQ(ierr); 3633b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 3643b224e63SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 3653b224e63SBarry Smith 3663b224e63SBarry Smith ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr); 3673b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr); 3683b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr); 3693b224e63SBarry Smith ilink = jac->head; 3703b224e63SBarry Smith ilink->x = jac->x[0]; ilink->y = jac->y[0]; 3713b224e63SBarry Smith ilink = ilink->next; 3723b224e63SBarry Smith ilink->x = jac->x[1]; ilink->y = jac->y[1]; 3733b224e63SBarry Smith } 3743b224e63SBarry Smith } else { 37597bbdb24SBarry Smith /* set up the individual PCs */ 37697bbdb24SBarry Smith i = 0; 3775a9f2f41SSatish Balay ilink = jac->head; 3785a9f2f41SSatish Balay while (ilink) { 3795a9f2f41SSatish Balay ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr); 3803b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 3815a9f2f41SSatish Balay ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr); 3825a9f2f41SSatish Balay ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 38397bbdb24SBarry Smith i++; 3845a9f2f41SSatish Balay ilink = ilink->next; 3850971522cSBarry Smith } 38697bbdb24SBarry Smith 38797bbdb24SBarry Smith /* create work vectors for each split */ 3881b9fc7fcSBarry Smith if (!jac->x) { 38997bbdb24SBarry Smith ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 3905a9f2f41SSatish Balay ilink = jac->head; 39197bbdb24SBarry Smith for (i=0; i<nsplit; i++) { 392906ed7ccSBarry Smith Vec *vl,*vr; 393906ed7ccSBarry Smith 394906ed7ccSBarry Smith ierr = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr); 395906ed7ccSBarry Smith ilink->x = *vr; 396906ed7ccSBarry Smith ilink->y = *vl; 397906ed7ccSBarry Smith ierr = PetscFree(vr);CHKERRQ(ierr); 398906ed7ccSBarry Smith ierr = PetscFree(vl);CHKERRQ(ierr); 3995a9f2f41SSatish Balay jac->x[i] = ilink->x; 4005a9f2f41SSatish Balay jac->y[i] = ilink->y; 4015a9f2f41SSatish Balay ilink = ilink->next; 40297bbdb24SBarry Smith } 4033b224e63SBarry Smith } 4043b224e63SBarry Smith } 4053b224e63SBarry Smith 4063b224e63SBarry Smith 407d0f46423SBarry Smith if (!jac->head->sctx) { 4083b224e63SBarry Smith Vec xtmp; 4093b224e63SBarry Smith 41079416396SBarry Smith /* compute scatter contexts needed by multiplicative versions and non-default splits */ 4111b9fc7fcSBarry Smith 4125a9f2f41SSatish Balay ilink = jac->head; 4131b9fc7fcSBarry Smith ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 4141b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 415704ba839SBarry Smith ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 4165a9f2f41SSatish Balay ilink = ilink->next; 4171b9fc7fcSBarry Smith } 4181b9fc7fcSBarry Smith ierr = VecDestroy(xtmp);CHKERRQ(ierr); 4191b9fc7fcSBarry Smith } 4200971522cSBarry Smith PetscFunctionReturn(0); 4210971522cSBarry Smith } 4220971522cSBarry Smith 4235a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 424ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 425ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 4265a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 427ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 428ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 42979416396SBarry Smith 4300971522cSBarry Smith #undef __FUNCT__ 4313b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 4323b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 4333b224e63SBarry Smith { 4343b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 4353b224e63SBarry Smith PetscErrorCode ierr; 4363b224e63SBarry Smith KSP ksp; 4373b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 4383b224e63SBarry Smith 4393b224e63SBarry Smith PetscFunctionBegin; 4403b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 4413b224e63SBarry Smith 4423b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4433b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4443b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 4453b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 4463b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 4473b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4483b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4493b224e63SBarry Smith 4503b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 4513b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4523b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4533b224e63SBarry Smith 4543b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 4553b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 4563b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 4573b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4583b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4593b224e63SBarry Smith 4603b224e63SBarry Smith PetscFunctionReturn(0); 4613b224e63SBarry Smith } 4623b224e63SBarry Smith 4633b224e63SBarry Smith #undef __FUNCT__ 4640971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 4650971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 4660971522cSBarry Smith { 4670971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 4680971522cSBarry Smith PetscErrorCode ierr; 4695a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 4703e197d65SBarry Smith PetscInt bs,cnt; 4710971522cSBarry Smith 4720971522cSBarry Smith PetscFunctionBegin; 47351f519a2SBarry Smith CHKMEMQ; 474e25d487eSBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 47551f519a2SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 47651f519a2SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 47751f519a2SBarry Smith 47879416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 4791b9fc7fcSBarry Smith if (jac->defaultsplit) { 4800971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 4815a9f2f41SSatish Balay while (ilink) { 4825a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 4835a9f2f41SSatish Balay ilink = ilink->next; 4840971522cSBarry Smith } 4850971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 4861b9fc7fcSBarry Smith } else { 487efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 4885a9f2f41SSatish Balay while (ilink) { 4895a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 4905a9f2f41SSatish Balay ilink = ilink->next; 4911b9fc7fcSBarry Smith } 4921b9fc7fcSBarry Smith } 49316913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 49479416396SBarry Smith if (!jac->w1) { 49579416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 49679416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 49779416396SBarry Smith } 498efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 4995a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 5003e197d65SBarry Smith cnt = 1; 5015a9f2f41SSatish Balay while (ilink->next) { 5025a9f2f41SSatish Balay ilink = ilink->next; 5033e197d65SBarry Smith if (jac->Afield) { 5043e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 5053e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 5063e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 5073e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5083e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5093e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 5103e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5113e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5123e197d65SBarry Smith } else { 5133e197d65SBarry Smith /* compute the residual over the entire vector */ 5141ab39975SBarry Smith ierr = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr); 515efb30889SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 5165a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 51779416396SBarry Smith } 5183e197d65SBarry Smith } 51951f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 52011755939SBarry Smith cnt -= 2; 52151f519a2SBarry Smith while (ilink->previous) { 52251f519a2SBarry Smith ilink = ilink->previous; 52311755939SBarry Smith if (jac->Afield) { 52411755939SBarry Smith /* compute the residual only over the part of the vector needed */ 52511755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 52611755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 52711755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 52811755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 52911755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 53011755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 53111755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 53211755939SBarry Smith } else { 5331ab39975SBarry Smith ierr = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr); 53451f519a2SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 53551f519a2SBarry Smith ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 53679416396SBarry Smith } 53751f519a2SBarry Smith } 53811755939SBarry Smith } 53916913363SBarry Smith } else SETERRQ1(PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 54051f519a2SBarry Smith CHKMEMQ; 5410971522cSBarry Smith PetscFunctionReturn(0); 5420971522cSBarry Smith } 5430971522cSBarry Smith 544421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 545ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 546ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 547421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 548ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 549ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 550421e10b8SBarry Smith 551421e10b8SBarry Smith #undef __FUNCT__ 552421e10b8SBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 553421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 554421e10b8SBarry Smith { 555421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 556421e10b8SBarry Smith PetscErrorCode ierr; 557421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 558421e10b8SBarry Smith PetscInt bs; 559421e10b8SBarry Smith 560421e10b8SBarry Smith PetscFunctionBegin; 561421e10b8SBarry Smith CHKMEMQ; 562421e10b8SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 563421e10b8SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 564421e10b8SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 565421e10b8SBarry Smith 566421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 567421e10b8SBarry Smith if (jac->defaultsplit) { 568421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 569421e10b8SBarry Smith while (ilink) { 570421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 571421e10b8SBarry Smith ilink = ilink->next; 572421e10b8SBarry Smith } 573421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 574421e10b8SBarry Smith } else { 575421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 576421e10b8SBarry Smith while (ilink) { 577421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 578421e10b8SBarry Smith ilink = ilink->next; 579421e10b8SBarry Smith } 580421e10b8SBarry Smith } 581421e10b8SBarry Smith } else { 582421e10b8SBarry Smith if (!jac->w1) { 583421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 584421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 585421e10b8SBarry Smith } 586421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 587421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 588421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 589421e10b8SBarry Smith while (ilink->next) { 590421e10b8SBarry Smith ilink = ilink->next; 5919989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 592421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 593421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 594421e10b8SBarry Smith } 595421e10b8SBarry Smith while (ilink->previous) { 596421e10b8SBarry Smith ilink = ilink->previous; 5979989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 598421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 599421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 600421e10b8SBarry Smith } 601421e10b8SBarry Smith } else { 602421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 603421e10b8SBarry Smith ilink = ilink->next; 604421e10b8SBarry Smith } 605421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 606421e10b8SBarry Smith while (ilink->previous) { 607421e10b8SBarry Smith ilink = ilink->previous; 6089989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 609421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 610421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 611421e10b8SBarry Smith } 612421e10b8SBarry Smith } 613421e10b8SBarry Smith } 614421e10b8SBarry Smith CHKMEMQ; 615421e10b8SBarry Smith PetscFunctionReturn(0); 616421e10b8SBarry Smith } 617421e10b8SBarry Smith 6180971522cSBarry Smith #undef __FUNCT__ 6190971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 6200971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 6210971522cSBarry Smith { 6220971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6230971522cSBarry Smith PetscErrorCode ierr; 6245a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 6250971522cSBarry Smith 6260971522cSBarry Smith PetscFunctionBegin; 6275a9f2f41SSatish Balay while (ilink) { 6285a9f2f41SSatish Balay ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr); 6295a9f2f41SSatish Balay if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);} 6305a9f2f41SSatish Balay if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);} 6315a9f2f41SSatish Balay if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);} 632704ba839SBarry Smith if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);} 6335a9f2f41SSatish Balay next = ilink->next; 634704ba839SBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 635704ba839SBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 6365a9f2f41SSatish Balay ilink = next; 6370971522cSBarry Smith } 63805b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 63997bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 64068dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 64179416396SBarry Smith if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);} 64279416396SBarry Smith if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);} 6433b224e63SBarry Smith if (jac->schur) {ierr = MatDestroy(jac->schur);CHKERRQ(ierr);} 644084e4875SJed Brown if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);} 645d0f46423SBarry Smith if (jac->kspschur) {ierr = KSPDestroy(jac->kspschur);CHKERRQ(ierr);} 6463b224e63SBarry Smith if (jac->B) {ierr = MatDestroy(jac->B);CHKERRQ(ierr);} 6473b224e63SBarry Smith if (jac->C) {ierr = MatDestroy(jac->C);CHKERRQ(ierr);} 6480971522cSBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 6490971522cSBarry Smith PetscFunctionReturn(0); 6500971522cSBarry Smith } 6510971522cSBarry Smith 6520971522cSBarry Smith #undef __FUNCT__ 6530971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 6540971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 6550971522cSBarry Smith { 6561b9fc7fcSBarry Smith PetscErrorCode ierr; 65751f519a2SBarry Smith PetscInt i = 0,nfields,*fields,bs; 658084e4875SJed Brown PetscTruth flg; 6591b9fc7fcSBarry Smith char optionname[128]; 6609dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6613b224e63SBarry Smith PCCompositeType ctype; 6621b9fc7fcSBarry Smith 6630971522cSBarry Smith PetscFunctionBegin; 6641b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 66551f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 66651f519a2SBarry Smith if (flg) { 66751f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 66851f519a2SBarry Smith } 669704ba839SBarry Smith 6703b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 6713b224e63SBarry Smith if (flg) { 6723b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 6733b224e63SBarry Smith } 674d32f9abdSBarry Smith 675c30613efSMatthew Knepley /* Only setup fields once */ 676c30613efSMatthew Knepley if ((jac->bs > 0) && (jac->nsplits == 0)) { 677d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 678d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 67951f519a2SBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr); 6801b9fc7fcSBarry Smith while (PETSC_TRUE) { 68113f74950SBarry Smith sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 68251f519a2SBarry Smith nfields = jac->bs; 6831b9fc7fcSBarry Smith ierr = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr); 6841b9fc7fcSBarry Smith if (!flg) break; 6851b9fc7fcSBarry Smith if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 6861b9fc7fcSBarry Smith ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr); 6871b9fc7fcSBarry Smith } 68851f519a2SBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 689d32f9abdSBarry Smith } 690084e4875SJed 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); 6911b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 6920971522cSBarry Smith PetscFunctionReturn(0); 6930971522cSBarry Smith } 6940971522cSBarry Smith 6950971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 6960971522cSBarry Smith 6970971522cSBarry Smith EXTERN_C_BEGIN 6980971522cSBarry Smith #undef __FUNCT__ 6990971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 700dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields) 7010971522cSBarry Smith { 70297bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7030971522cSBarry Smith PetscErrorCode ierr; 7045a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 70569a612a9SBarry Smith char prefix[128]; 70651f519a2SBarry Smith PetscInt i; 7070971522cSBarry Smith 7080971522cSBarry Smith PetscFunctionBegin; 7090971522cSBarry Smith if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested"); 71051f519a2SBarry Smith for (i=0; i<n; i++) { 71151f519a2SBarry Smith if (fields[i] >= jac->bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs); 71251f519a2SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 71351f519a2SBarry Smith } 714704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 715704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 7165a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 7175a9f2f41SSatish Balay ilink->nfields = n; 7185a9f2f41SSatish Balay ilink->next = PETSC_NULL; 7197adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 7201cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 7215a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 72269a612a9SBarry Smith 7237adad957SLisandro Dalcin if (((PetscObject)pc)->prefix) { 7247adad957SLisandro Dalcin sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 72569a612a9SBarry Smith } else { 72613f74950SBarry Smith sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 72769a612a9SBarry Smith } 7285a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 7290971522cSBarry Smith 7300971522cSBarry Smith if (!next) { 7315a9f2f41SSatish Balay jac->head = ilink; 73251f519a2SBarry Smith ilink->previous = PETSC_NULL; 7330971522cSBarry Smith } else { 7340971522cSBarry Smith while (next->next) { 7350971522cSBarry Smith next = next->next; 7360971522cSBarry Smith } 7375a9f2f41SSatish Balay next->next = ilink; 73851f519a2SBarry Smith ilink->previous = next; 7390971522cSBarry Smith } 7400971522cSBarry Smith jac->nsplits++; 7410971522cSBarry Smith PetscFunctionReturn(0); 7420971522cSBarry Smith } 7430971522cSBarry Smith EXTERN_C_END 7440971522cSBarry Smith 745e69d4d44SBarry Smith EXTERN_C_BEGIN 746e69d4d44SBarry Smith #undef __FUNCT__ 747e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 748e69d4d44SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 749e69d4d44SBarry Smith { 750e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 751e69d4d44SBarry Smith PetscErrorCode ierr; 752e69d4d44SBarry Smith 753e69d4d44SBarry Smith PetscFunctionBegin; 754e69d4d44SBarry Smith ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 755e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 756e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 757084e4875SJed Brown *n = jac->nsplits; 758e69d4d44SBarry Smith PetscFunctionReturn(0); 759e69d4d44SBarry Smith } 760e69d4d44SBarry Smith EXTERN_C_END 7610971522cSBarry Smith 7620971522cSBarry Smith EXTERN_C_BEGIN 7630971522cSBarry Smith #undef __FUNCT__ 76469a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 765dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 7660971522cSBarry Smith { 7670971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7680971522cSBarry Smith PetscErrorCode ierr; 7690971522cSBarry Smith PetscInt cnt = 0; 7705a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 7710971522cSBarry Smith 7720971522cSBarry Smith PetscFunctionBegin; 77369a612a9SBarry Smith ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 7745a9f2f41SSatish Balay while (ilink) { 7755a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 7765a9f2f41SSatish Balay ilink = ilink->next; 7770971522cSBarry Smith } 7780971522cSBarry Smith if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits); 7790971522cSBarry Smith *n = jac->nsplits; 7800971522cSBarry Smith PetscFunctionReturn(0); 7810971522cSBarry Smith } 7820971522cSBarry Smith EXTERN_C_END 7830971522cSBarry Smith 784704ba839SBarry Smith EXTERN_C_BEGIN 785704ba839SBarry Smith #undef __FUNCT__ 786704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 787704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,IS is) 788704ba839SBarry Smith { 789704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 790704ba839SBarry Smith PetscErrorCode ierr; 791704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 792704ba839SBarry Smith char prefix[128]; 793704ba839SBarry Smith 794704ba839SBarry Smith PetscFunctionBegin; 79516913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 7961ab39975SBarry Smith ilink->is = is; 7971ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 798704ba839SBarry Smith ilink->next = PETSC_NULL; 799704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 8001cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 801704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 802704ba839SBarry Smith 803704ba839SBarry Smith if (((PetscObject)pc)->prefix) { 804704ba839SBarry Smith sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 805704ba839SBarry Smith } else { 806704ba839SBarry Smith sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 807704ba839SBarry Smith } 808704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 809704ba839SBarry Smith 810704ba839SBarry Smith if (!next) { 811704ba839SBarry Smith jac->head = ilink; 812704ba839SBarry Smith ilink->previous = PETSC_NULL; 813704ba839SBarry Smith } else { 814704ba839SBarry Smith while (next->next) { 815704ba839SBarry Smith next = next->next; 816704ba839SBarry Smith } 817704ba839SBarry Smith next->next = ilink; 818704ba839SBarry Smith ilink->previous = next; 819704ba839SBarry Smith } 820704ba839SBarry Smith jac->nsplits++; 821704ba839SBarry Smith 822704ba839SBarry Smith PetscFunctionReturn(0); 823704ba839SBarry Smith } 824704ba839SBarry Smith EXTERN_C_END 825704ba839SBarry Smith 8260971522cSBarry Smith #undef __FUNCT__ 8270971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 8280971522cSBarry Smith /*@ 8290971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 8300971522cSBarry Smith 8310971522cSBarry Smith Collective on PC 8320971522cSBarry Smith 8330971522cSBarry Smith Input Parameters: 8340971522cSBarry Smith + pc - the preconditioner context 8350971522cSBarry Smith . n - the number of fields in this split 8360971522cSBarry Smith . fields - the fields in this split 8370971522cSBarry Smith 8380971522cSBarry Smith Level: intermediate 8390971522cSBarry Smith 840d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 841d32f9abdSBarry Smith 842d32f9abdSBarry Smith The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block 843d32f9abdSBarry 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 844d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 845d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 846d32f9abdSBarry Smith 847d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 8480971522cSBarry Smith 8490971522cSBarry Smith @*/ 850dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields) 8510971522cSBarry Smith { 8520971522cSBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *); 8530971522cSBarry Smith 8540971522cSBarry Smith PetscFunctionBegin; 8550971522cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 8560971522cSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr); 8570971522cSBarry Smith if (f) { 8580971522cSBarry Smith ierr = (*f)(pc,n,fields);CHKERRQ(ierr); 8590971522cSBarry Smith } 8600971522cSBarry Smith PetscFunctionReturn(0); 8610971522cSBarry Smith } 8620971522cSBarry Smith 8630971522cSBarry Smith #undef __FUNCT__ 864704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 865704ba839SBarry Smith /*@ 866704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 867704ba839SBarry Smith 868704ba839SBarry Smith Collective on PC 869704ba839SBarry Smith 870704ba839SBarry Smith Input Parameters: 871704ba839SBarry Smith + pc - the preconditioner context 872704ba839SBarry Smith . is - the index set that defines the vector elements in this field 873704ba839SBarry Smith 874d32f9abdSBarry Smith 875d32f9abdSBarry Smith Notes: Use PCFieldSplitSetFields(), for fields defined by strided types. 876d32f9abdSBarry Smith 877704ba839SBarry Smith Level: intermediate 878704ba839SBarry Smith 879704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 880704ba839SBarry Smith 881704ba839SBarry Smith @*/ 882704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,IS is) 883704ba839SBarry Smith { 884704ba839SBarry Smith PetscErrorCode ierr,(*f)(PC,IS); 885704ba839SBarry Smith 886704ba839SBarry Smith PetscFunctionBegin; 887704ba839SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 888704ba839SBarry Smith PetscValidHeaderSpecific(is,IS_COOKIE,1); 889704ba839SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr); 890704ba839SBarry Smith if (f) { 891704ba839SBarry Smith ierr = (*f)(pc,is);CHKERRQ(ierr); 892704ba839SBarry Smith } 893704ba839SBarry Smith PetscFunctionReturn(0); 894704ba839SBarry Smith } 895704ba839SBarry Smith 896704ba839SBarry Smith #undef __FUNCT__ 89751f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 89851f519a2SBarry Smith /*@ 89951f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 90051f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 90151f519a2SBarry Smith 90251f519a2SBarry Smith Collective on PC 90351f519a2SBarry Smith 90451f519a2SBarry Smith Input Parameters: 90551f519a2SBarry Smith + pc - the preconditioner context 90651f519a2SBarry Smith - bs - the block size 90751f519a2SBarry Smith 90851f519a2SBarry Smith Level: intermediate 90951f519a2SBarry Smith 91051f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 91151f519a2SBarry Smith 91251f519a2SBarry Smith @*/ 91351f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 91451f519a2SBarry Smith { 91551f519a2SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt); 91651f519a2SBarry Smith 91751f519a2SBarry Smith PetscFunctionBegin; 91851f519a2SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 91951f519a2SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr); 92051f519a2SBarry Smith if (f) { 92151f519a2SBarry Smith ierr = (*f)(pc,bs);CHKERRQ(ierr); 92251f519a2SBarry Smith } 92351f519a2SBarry Smith PetscFunctionReturn(0); 92451f519a2SBarry Smith } 92551f519a2SBarry Smith 92651f519a2SBarry Smith #undef __FUNCT__ 92769a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 9280971522cSBarry Smith /*@C 92969a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 9300971522cSBarry Smith 93169a612a9SBarry Smith Collective on KSP 9320971522cSBarry Smith 9330971522cSBarry Smith Input Parameter: 9340971522cSBarry Smith . pc - the preconditioner context 9350971522cSBarry Smith 9360971522cSBarry Smith Output Parameters: 9370971522cSBarry Smith + n - the number of split 93869a612a9SBarry Smith - pc - the array of KSP contexts 9390971522cSBarry Smith 9400971522cSBarry Smith Note: 941d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 942d32f9abdSBarry Smith (not the KSP just the array that contains them). 9430971522cSBarry Smith 94469a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 9450971522cSBarry Smith 9460971522cSBarry Smith Level: advanced 9470971522cSBarry Smith 9480971522cSBarry Smith .seealso: PCFIELDSPLIT 9490971522cSBarry Smith @*/ 950dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 9510971522cSBarry Smith { 95269a612a9SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **); 9530971522cSBarry Smith 9540971522cSBarry Smith PetscFunctionBegin; 9550971522cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 9560971522cSBarry Smith PetscValidIntPointer(n,2); 95769a612a9SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 9580971522cSBarry Smith if (f) { 95969a612a9SBarry Smith ierr = (*f)(pc,n,subksp);CHKERRQ(ierr); 9600971522cSBarry Smith } else { 96169a612a9SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC"); 9620971522cSBarry Smith } 9630971522cSBarry Smith PetscFunctionReturn(0); 9640971522cSBarry Smith } 9650971522cSBarry Smith 966e69d4d44SBarry Smith #undef __FUNCT__ 967e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 968e69d4d44SBarry Smith /*@ 969e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 970e69d4d44SBarry Smith D matrix. Otherwise no preconditioner is used. 971e69d4d44SBarry Smith 972e69d4d44SBarry Smith Collective on PC 973e69d4d44SBarry Smith 974e69d4d44SBarry Smith Input Parameters: 975e69d4d44SBarry Smith + pc - the preconditioner context 976084e4875SJed Brown . ptype - which matrix to use for preconditioning the Schur complement 977084e4875SJed Brown - userpre - matrix to use for preconditioning, or PETSC_NULL 978084e4875SJed Brown 979084e4875SJed Brown Notes: 980084e4875SJed Brown The default is to use the block on the diagonal of the preconditioning matrix. This is D, in the (1,1) position. 981084e4875SJed Brown There are currently no preconditioners that work directly with the Schur complement so setting 982084e4875SJed Brown PC_FIELDSPLIT_SCHUR_PRE_SELF is observationally equivalent to -fieldsplit_1_pc_type none. 983e69d4d44SBarry Smith 984e69d4d44SBarry Smith Options Database: 985084e4875SJed Brown . -pc_fieldsplit_schur_precondition <self,user,diag> default is diag 986e69d4d44SBarry Smith 987e69d4d44SBarry Smith Level: intermediate 988e69d4d44SBarry Smith 989084e4875SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType 990e69d4d44SBarry Smith 991e69d4d44SBarry Smith @*/ 992084e4875SJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 993e69d4d44SBarry Smith { 994084e4875SJed Brown PetscErrorCode ierr,(*f)(PC,PCFieldSplitSchurPreType,Mat); 995e69d4d44SBarry Smith 996e69d4d44SBarry Smith PetscFunctionBegin; 997e69d4d44SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 998e69d4d44SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",(void (**)(void))&f);CHKERRQ(ierr); 999e69d4d44SBarry Smith if (f) { 1000084e4875SJed Brown ierr = (*f)(pc,ptype,pre);CHKERRQ(ierr); 1001e69d4d44SBarry Smith } 1002e69d4d44SBarry Smith PetscFunctionReturn(0); 1003e69d4d44SBarry Smith } 1004e69d4d44SBarry Smith 1005e69d4d44SBarry Smith EXTERN_C_BEGIN 1006e69d4d44SBarry Smith #undef __FUNCT__ 1007e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 1008084e4875SJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1009e69d4d44SBarry Smith { 1010e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1011084e4875SJed Brown PetscErrorCode ierr; 1012e69d4d44SBarry Smith 1013e69d4d44SBarry Smith PetscFunctionBegin; 1014084e4875SJed Brown jac->schurpre = ptype; 1015084e4875SJed Brown if (pre) { 1016084e4875SJed Brown if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);} 1017084e4875SJed Brown jac->schur_user = pre; 1018084e4875SJed Brown ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1019084e4875SJed Brown } 1020e69d4d44SBarry Smith PetscFunctionReturn(0); 1021e69d4d44SBarry Smith } 1022e69d4d44SBarry Smith EXTERN_C_END 1023e69d4d44SBarry Smith 102430ad9308SMatthew Knepley #undef __FUNCT__ 102530ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 102630ad9308SMatthew Knepley /*@C 102730ad9308SMatthew Knepley PCFieldSplitGetSchurBlocks - Gets the all matrix blocks for the Schur complement 102830ad9308SMatthew Knepley 102930ad9308SMatthew Knepley Collective on KSP 103030ad9308SMatthew Knepley 103130ad9308SMatthew Knepley Input Parameter: 103230ad9308SMatthew Knepley . pc - the preconditioner context 103330ad9308SMatthew Knepley 103430ad9308SMatthew Knepley Output Parameters: 103530ad9308SMatthew Knepley + A - the (0,0) block 103630ad9308SMatthew Knepley . B - the (0,1) block 103730ad9308SMatthew Knepley . C - the (1,0) block 103830ad9308SMatthew Knepley - D - the (1,1) block 103930ad9308SMatthew Knepley 104030ad9308SMatthew Knepley Level: advanced 104130ad9308SMatthew Knepley 104230ad9308SMatthew Knepley .seealso: PCFIELDSPLIT 104330ad9308SMatthew Knepley @*/ 104430ad9308SMatthew Knepley PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSchurBlocks(PC pc,Mat *A,Mat *B,Mat *C, Mat *D) 104530ad9308SMatthew Knepley { 104630ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 104730ad9308SMatthew Knepley 104830ad9308SMatthew Knepley PetscFunctionBegin; 104930ad9308SMatthew Knepley PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1050cf8ad460SMatthew Knepley if (jac->type != PC_COMPOSITE_SCHUR) {SETERRQ(PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");} 105130ad9308SMatthew Knepley if (A) *A = jac->pmat[0]; 105230ad9308SMatthew Knepley if (B) *B = jac->B; 105330ad9308SMatthew Knepley if (C) *C = jac->C; 105430ad9308SMatthew Knepley if (D) *D = jac->pmat[1]; 105530ad9308SMatthew Knepley PetscFunctionReturn(0); 105630ad9308SMatthew Knepley } 105730ad9308SMatthew Knepley 105879416396SBarry Smith EXTERN_C_BEGIN 105979416396SBarry Smith #undef __FUNCT__ 106079416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 1061dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 106279416396SBarry Smith { 106379416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1064e69d4d44SBarry Smith PetscErrorCode ierr; 106579416396SBarry Smith 106679416396SBarry Smith PetscFunctionBegin; 106779416396SBarry Smith jac->type = type; 10683b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 10693b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 10703b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 1071e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1072e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1073e69d4d44SBarry Smith 10743b224e63SBarry Smith } else { 10753b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 10763b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 1077e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 10789e7d6b0aSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 10793b224e63SBarry Smith } 108079416396SBarry Smith PetscFunctionReturn(0); 108179416396SBarry Smith } 108279416396SBarry Smith EXTERN_C_END 108379416396SBarry Smith 108451f519a2SBarry Smith EXTERN_C_BEGIN 108551f519a2SBarry Smith #undef __FUNCT__ 108651f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 108751f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 108851f519a2SBarry Smith { 108951f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 109051f519a2SBarry Smith 109151f519a2SBarry Smith PetscFunctionBegin; 109251f519a2SBarry Smith if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 109351f519a2SBarry Smith if (jac->bs > 0 && jac->bs != bs) SETERRQ2(PETSC_ERR_ARG_WRONGSTATE,"Cannot change fieldsplit blocksize from %D to %D after it has been set",jac->bs,bs); 109451f519a2SBarry Smith jac->bs = bs; 109551f519a2SBarry Smith PetscFunctionReturn(0); 109651f519a2SBarry Smith } 109751f519a2SBarry Smith EXTERN_C_END 109851f519a2SBarry Smith 109979416396SBarry Smith #undef __FUNCT__ 110079416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1101bc08b0f1SBarry Smith /*@ 110279416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 110379416396SBarry Smith 110479416396SBarry Smith Collective on PC 110579416396SBarry Smith 110679416396SBarry Smith Input Parameter: 110779416396SBarry Smith . pc - the preconditioner context 110881540f2fSBarry Smith . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 110979416396SBarry Smith 111079416396SBarry Smith Options Database Key: 1111a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 111279416396SBarry Smith 111379416396SBarry Smith Level: Developer 111479416396SBarry Smith 111579416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 111679416396SBarry Smith 111779416396SBarry Smith .seealso: PCCompositeSetType() 111879416396SBarry Smith 111979416396SBarry Smith @*/ 1120dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type) 112179416396SBarry Smith { 112279416396SBarry Smith PetscErrorCode ierr,(*f)(PC,PCCompositeType); 112379416396SBarry Smith 112479416396SBarry Smith PetscFunctionBegin; 112579416396SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 112679416396SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 112779416396SBarry Smith if (f) { 112879416396SBarry Smith ierr = (*f)(pc,type);CHKERRQ(ierr); 112979416396SBarry Smith } 113079416396SBarry Smith PetscFunctionReturn(0); 113179416396SBarry Smith } 113279416396SBarry Smith 11330971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 11340971522cSBarry Smith /*MC 1135a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 11360971522cSBarry Smith fields or groups of fields 11370971522cSBarry Smith 11380971522cSBarry Smith 1139edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1140edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 11410971522cSBarry Smith 1142a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 114369a612a9SBarry Smith and set the options directly on the resulting KSP object 11440971522cSBarry Smith 11450971522cSBarry Smith Level: intermediate 11460971522cSBarry Smith 114779416396SBarry Smith Options Database Keys: 114881540f2fSBarry Smith + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 114981540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 115081540f2fSBarry Smith been supplied explicitly by -pc_fieldsplit_%d_fields 115181540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 115281540f2fSBarry Smith . -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative> 1153e69d4d44SBarry Smith . -pc_fieldsplit_schur_precondition <true,false> default is true 115479416396SBarry Smith 1155edf189efSBarry Smith - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 11563b224e63SBarry Smith for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 11573b224e63SBarry Smith 11583b224e63SBarry Smith 1159d32f9abdSBarry Smith Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1160d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1161d32f9abdSBarry Smith 1162d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1163d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1164d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1165d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1166d32f9abdSBarry Smith 1167d32f9abdSBarry Smith Currently for the multiplicative version, the updated residual needed for the next field 1168d32f9abdSBarry Smith solve is computed via a matrix vector product over the entire array. An optimization would be 1169d32f9abdSBarry Smith to update the residual only for the part of the right hand side associated with the next field 1170d32f9abdSBarry Smith solve. (This would involve more MatGetSubMatrix() calls or some other mechanism to compute the 1171d32f9abdSBarry Smith part of the matrix needed to just update part of the residual). 1172d32f9abdSBarry Smith 1173e69d4d44SBarry Smith For the Schur complement preconditioner if J = ( A B ) 1174e69d4d44SBarry Smith ( C D ) 1175e69d4d44SBarry Smith the preconditioner is 1176e69d4d44SBarry Smith (I -B inv(A)) ( inv(A) 0 ) (I 0 ) 1177e69d4d44SBarry Smith (0 I ) ( 0 inv(S) ) (-C inv(A) I ) 1178edf189efSBarry Smith where the action of inv(A) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1179e69d4d44SBarry Smith inv(S) is computed using the KSP solver with prefix -schur_. For PCFieldSplitGetKSP() when field number is 1180edf189efSBarry Smith 0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1181e69d4d44SBarry Smith D is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1182e69d4d44SBarry Smith option. 1183e69d4d44SBarry Smith 1184edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1185edf189efSBarry Smith is used automatically for a second block. 1186edf189efSBarry Smith 1187a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 11880971522cSBarry Smith 1189a541d17aSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners 1190e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 11910971522cSBarry Smith M*/ 11920971522cSBarry Smith 11930971522cSBarry Smith EXTERN_C_BEGIN 11940971522cSBarry Smith #undef __FUNCT__ 11950971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 1196dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc) 11970971522cSBarry Smith { 11980971522cSBarry Smith PetscErrorCode ierr; 11990971522cSBarry Smith PC_FieldSplit *jac; 12000971522cSBarry Smith 12010971522cSBarry Smith PetscFunctionBegin; 120238f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 12030971522cSBarry Smith jac->bs = -1; 12040971522cSBarry Smith jac->nsplits = 0; 12053e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1206*e6cab6aaSJed Brown jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 120751f519a2SBarry Smith 12080971522cSBarry Smith pc->data = (void*)jac; 12090971522cSBarry Smith 12100971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1211421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 12120971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 12130971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 12140971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 12150971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 12160971522cSBarry Smith pc->ops->applyrichardson = 0; 12170971522cSBarry Smith 121869a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 121969a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 12200971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 12210971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1222704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1223704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 122479416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 122579416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 122651f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 122751f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 12280971522cSBarry Smith PetscFunctionReturn(0); 12290971522cSBarry Smith } 12300971522cSBarry Smith EXTERN_C_END 12310971522cSBarry Smith 12320971522cSBarry Smith 1233a541d17aSBarry Smith /*MC 1234a541d17aSBarry Smith Block_Preconditioners - PETSc provides support for a variety of "block preconditioners", this provides an 1235a541d17aSBarry Smith overview of these methods. 1236a541d17aSBarry Smith 1237a541d17aSBarry Smith Consider the solution to ( A B ) (x_1) = (b_1) 1238a541d17aSBarry Smith ( C D ) (x_2) (b_2) 1239a541d17aSBarry Smith 1240a541d17aSBarry Smith Important special cases, the Stokes equation: C = B' and D = 0 (A B) (x_1) = (b_1) 1241a541d17aSBarry Smith B' 0) (x_2) (b_2) 1242a541d17aSBarry Smith 1243a541d17aSBarry Smith One of the goals of the PCFieldSplit preconditioner in PETSc is to provide a variety of preconditioners 1244a541d17aSBarry Smith for this block system. 1245a541d17aSBarry Smith 1246a541d17aSBarry Smith Consider an additional matrix (Ap Bp) 1247a541d17aSBarry Smith (Cp Dp) where some or all of the entries may be the same as 1248a541d17aSBarry Smith in the original matrix (for example Ap == A). 1249a541d17aSBarry Smith 1250a541d17aSBarry Smith In the following, A^ denotes the approximate application of the inverse of A, possibly using Ap in the 1251a541d17aSBarry Smith approximation. In PETSc this simply means one has called KSPSetOperators(ksp,A,Ap,...) or KSPSetOperators(ksp,Ap,Ap,...) 1252a541d17aSBarry Smith 1253a541d17aSBarry Smith Block Jacobi: x_1 = A^ b_1 1254a541d17aSBarry Smith x_2 = D^ b_2 1255a541d17aSBarry Smith 1256a541d17aSBarry Smith Lower block Gauss-Seidel: x_1 = A^ b_1 1257a541d17aSBarry Smith x_2 = D^ (b_2 - C x_1) variant x_2 = D^ (b_2 - Cp x_1) 1258a541d17aSBarry Smith 1259a541d17aSBarry 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) 1260a541d17aSBarry Smith Interestingly this form is not actually a symmetric matrix, the symmetric version is 1261a541d17aSBarry Smith x_1 = A^(b_1 - B x_2) variant x_1 = A^(b_1 - Bp x_2) 1262a541d17aSBarry Smith 12630bc0a719SBarry Smith Level: intermediate 12640bc0a719SBarry Smith 1265a541d17aSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCFIELDSPLIT 1266a541d17aSBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 1267a541d17aSBarry Smith M*/ 1268