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 8*084e4875SJed Brown const char *PCFieldSplitSchurPreTypes[] = {"self","diag","user","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0}; 9*084e4875SJed 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; 17704ba839SBarry Smith IS is,cis; 18704ba839SBarry Smith PetscInt csize; 1951f519a2SBarry Smith PC_FieldSplitLink next,previous; 200971522cSBarry Smith }; 210971522cSBarry Smith 220971522cSBarry Smith typedef struct { 2368dd23aaSBarry Smith PCCompositeType type; 24a4efd8eaSMatthew Knepley PetscTruth defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */ 2530ad9308SMatthew Knepley PetscInt bs; /* Block size for IS and Mat structures */ 2630ad9308SMatthew Knepley PetscInt nsplits; /* Number of field divisions defined */ 2779416396SBarry Smith Vec *x,*y,w1,w2; 2830ad9308SMatthew Knepley Mat *pmat; /* The diagonal block for each split */ 2930ad9308SMatthew Knepley Mat *Afield; /* The rows of the matrix associated with each split */ 30704ba839SBarry Smith PetscTruth issetup; 3130ad9308SMatthew Knepley /* Only used when Schur complement preconditioning is used */ 3230ad9308SMatthew Knepley Mat B; /* The (0,1) block */ 3330ad9308SMatthew Knepley Mat C; /* The (1,0) block */ 3430ad9308SMatthew Knepley Mat schur; /* The Schur complement S = D - C A^{-1} B */ 35*084e4875SJed Brown Mat schur_user; /* User-provided preconditioning matrix for the Schur complement */ 36*084e4875SJed Brown PCFieldSplitSchurPreType schurpre; /* Determines which preconditioning matrix is used for the Schur complement */ 3730ad9308SMatthew Knepley KSP kspschur; /* The solver for S */ 3897bbdb24SBarry Smith PC_FieldSplitLink head; 390971522cSBarry Smith } PC_FieldSplit; 400971522cSBarry Smith 4116913363SBarry Smith /* 4216913363SBarry Smith Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of 4316913363SBarry Smith inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the 4416913363SBarry Smith PC you could change this. 4516913363SBarry Smith */ 46*084e4875SJed Brown 47*084e4875SJed Brown /* This helper is so that setting a user-provided preconditioning matrix orthogonal to choosing to use it. This way the 48*084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */ 49*084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac) 50*084e4875SJed Brown { 51*084e4875SJed Brown switch (jac->schurpre) { 52*084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur; 53*084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1]; 54*084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */ 55*084e4875SJed Brown default: 56*084e4875SJed Brown return jac->schur_user ? jac->schur_user : jac->pmat[1]; 57*084e4875SJed Brown } 58*084e4875SJed Brown } 59*084e4875SJed Brown 60*084e4875SJed Brown 610971522cSBarry Smith #undef __FUNCT__ 620971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit" 630971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 640971522cSBarry Smith { 650971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 660971522cSBarry Smith PetscErrorCode ierr; 670971522cSBarry Smith PetscTruth iascii; 680971522cSBarry Smith PetscInt i,j; 695a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 700971522cSBarry Smith 710971522cSBarry Smith PetscFunctionBegin; 720971522cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 730971522cSBarry Smith if (iascii) { 7451f519a2SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr); 7569a612a9SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 760971522cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 770971522cSBarry Smith for (i=0; i<jac->nsplits; i++) { 781ab39975SBarry Smith if (ilink->fields) { 790971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 8079416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 815a9f2f41SSatish Balay for (j=0; j<ilink->nfields; j++) { 8279416396SBarry Smith if (j > 0) { 8379416396SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 8479416396SBarry Smith } 855a9f2f41SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 860971522cSBarry Smith } 870971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 8879416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 891ab39975SBarry Smith } else { 901ab39975SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 911ab39975SBarry Smith } 925a9f2f41SSatish Balay ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 935a9f2f41SSatish Balay ilink = ilink->next; 940971522cSBarry Smith } 950971522cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 960971522cSBarry Smith } else { 970971522cSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 980971522cSBarry Smith } 990971522cSBarry Smith PetscFunctionReturn(0); 1000971522cSBarry Smith } 1010971522cSBarry Smith 1020971522cSBarry Smith #undef __FUNCT__ 1033b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur" 1043b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer) 1053b224e63SBarry Smith { 1063b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1073b224e63SBarry Smith PetscErrorCode ierr; 1083b224e63SBarry Smith PetscTruth iascii; 1093b224e63SBarry Smith PetscInt i,j; 1103b224e63SBarry Smith PC_FieldSplitLink ilink = jac->head; 1113b224e63SBarry Smith KSP ksp; 1123b224e63SBarry Smith 1133b224e63SBarry Smith PetscFunctionBegin; 1143b224e63SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1153b224e63SBarry Smith if (iascii) { 1163b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, blocksize = %D\n",jac->bs);CHKERRQ(ierr); 1173b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Split info:\n");CHKERRQ(ierr); 1183b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1193b224e63SBarry Smith for (i=0; i<jac->nsplits; i++) { 1203b224e63SBarry Smith if (ilink->fields) { 1213b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 1223b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 1233b224e63SBarry Smith for (j=0; j<ilink->nfields; j++) { 1243b224e63SBarry Smith if (j > 0) { 1253b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 1263b224e63SBarry Smith } 1273b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 1283b224e63SBarry Smith } 1293b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1303b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1313b224e63SBarry Smith } else { 1323b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1333b224e63SBarry Smith } 1343b224e63SBarry Smith ilink = ilink->next; 1353b224e63SBarry Smith } 1363b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A block \n");CHKERRQ(ierr); 1373b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 1383b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1393b224e63SBarry Smith ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 1403b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1413b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = D - C inv(A) B \n");CHKERRQ(ierr); 1423b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1433b224e63SBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 1443b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1453b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1463b224e63SBarry Smith } else { 1473b224e63SBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1483b224e63SBarry Smith } 1493b224e63SBarry Smith PetscFunctionReturn(0); 1503b224e63SBarry Smith } 1513b224e63SBarry Smith 1523b224e63SBarry Smith #undef __FUNCT__ 15369a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 15469a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 1550971522cSBarry Smith { 1560971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1570971522cSBarry Smith PetscErrorCode ierr; 1585a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 159d32f9abdSBarry Smith PetscInt i = 0,*ifields,nfields; 160d32f9abdSBarry Smith PetscTruth flg = PETSC_FALSE,*fields,flg2; 161d32f9abdSBarry Smith char optionname[128]; 1620971522cSBarry Smith 1630971522cSBarry Smith PetscFunctionBegin; 164d32f9abdSBarry Smith if (!ilink) { 165d32f9abdSBarry Smith 166521d7252SBarry Smith if (jac->bs <= 0) { 167704ba839SBarry Smith if (pc->pmat) { 168521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 169704ba839SBarry Smith } else { 170704ba839SBarry Smith jac->bs = 1; 171704ba839SBarry Smith } 172521d7252SBarry Smith } 173d32f9abdSBarry Smith 174d32f9abdSBarry Smith ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 175d32f9abdSBarry Smith if (!flg) { 176d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 177d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 178d32f9abdSBarry Smith flg = PETSC_TRUE; /* switched off automatically if user sets fields manually here */ 179d32f9abdSBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 180d32f9abdSBarry Smith while (PETSC_TRUE) { 181d32f9abdSBarry Smith sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 182d32f9abdSBarry Smith nfields = jac->bs; 183d32f9abdSBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg2);CHKERRQ(ierr); 184d32f9abdSBarry Smith if (!flg2) break; 185d32f9abdSBarry Smith if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 186d32f9abdSBarry Smith flg = PETSC_FALSE; 187d32f9abdSBarry Smith ierr = PCFieldSplitSetFields(pc,nfields,ifields);CHKERRQ(ierr); 188d32f9abdSBarry Smith } 189d32f9abdSBarry Smith ierr = PetscFree(ifields);CHKERRQ(ierr); 190d32f9abdSBarry Smith } 191d32f9abdSBarry Smith 192d32f9abdSBarry Smith if (flg) { 193d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 19479416396SBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr); 19579416396SBarry Smith ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr); 1965a9f2f41SSatish Balay while (ilink) { 1975a9f2f41SSatish Balay for (i=0; i<ilink->nfields; i++) { 1985a9f2f41SSatish Balay fields[ilink->fields[i]] = PETSC_TRUE; 199521d7252SBarry Smith } 2005a9f2f41SSatish Balay ilink = ilink->next; 20179416396SBarry Smith } 20297bbdb24SBarry Smith jac->defaultsplit = PETSC_TRUE; 20379416396SBarry Smith for (i=0; i<jac->bs; i++) { 20479416396SBarry Smith if (!fields[i]) { 20579416396SBarry Smith ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr); 20679416396SBarry Smith } else { 20779416396SBarry Smith jac->defaultsplit = PETSC_FALSE; 20879416396SBarry Smith } 20979416396SBarry Smith } 21068dd23aaSBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 211521d7252SBarry Smith } 212edf189efSBarry Smith } else if (jac->nsplits == 1) { 213edf189efSBarry Smith if (ilink->is) { 214edf189efSBarry Smith IS is2; 215edf189efSBarry Smith PetscInt nmin,nmax; 216edf189efSBarry Smith 217edf189efSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 218edf189efSBarry Smith ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 219edf189efSBarry Smith ierr = PCFieldSplitSetIS(pc,is2);CHKERRQ(ierr); 220edf189efSBarry Smith ierr = ISDestroy(is2);CHKERRQ(ierr); 221edf189efSBarry Smith } else { 222edf189efSBarry Smith SETERRQ(PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 223edf189efSBarry Smith } 224d32f9abdSBarry Smith } 22569a612a9SBarry Smith PetscFunctionReturn(0); 22669a612a9SBarry Smith } 22769a612a9SBarry Smith 22869a612a9SBarry Smith 22969a612a9SBarry Smith #undef __FUNCT__ 23069a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 23169a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 23269a612a9SBarry Smith { 23369a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 23469a612a9SBarry Smith PetscErrorCode ierr; 2355a9f2f41SSatish Balay PC_FieldSplitLink ilink; 236cf502942SBarry Smith PetscInt i,nsplit,ccsize; 23769a612a9SBarry Smith MatStructure flag = pc->flag; 23868dd23aaSBarry Smith PetscTruth sorted,getsub; 23969a612a9SBarry Smith 24069a612a9SBarry Smith PetscFunctionBegin; 24169a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 24297bbdb24SBarry Smith nsplit = jac->nsplits; 2435a9f2f41SSatish Balay ilink = jac->head; 24497bbdb24SBarry Smith 24597bbdb24SBarry Smith /* get the matrices for each split */ 246704ba839SBarry Smith if (!jac->issetup) { 2471b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 24897bbdb24SBarry Smith 249704ba839SBarry Smith jac->issetup = PETSC_TRUE; 250704ba839SBarry Smith 251704ba839SBarry Smith /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 25251f519a2SBarry Smith bs = jac->bs; 25397bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 254cf502942SBarry Smith ierr = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr); 2551b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 2561b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 2571b9fc7fcSBarry Smith if (jac->defaultsplit) { 258704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 259704ba839SBarry Smith } else if (!ilink->is) { 260ccb205f8SBarry Smith if (ilink->nfields > 1) { 2615a9f2f41SSatish Balay PetscInt *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields; 2625a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 2631b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 2641b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 2651b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 26697bbdb24SBarry Smith } 26797bbdb24SBarry Smith } 268704ba839SBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr); 269ccb205f8SBarry Smith ierr = PetscFree(ii);CHKERRQ(ierr); 270ccb205f8SBarry Smith } else { 271704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 272ccb205f8SBarry Smith } 2733e197d65SBarry Smith } 2743e197d65SBarry Smith ierr = ISGetLocalSize(ilink->is,&ilink->csize);CHKERRQ(ierr); 275704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 2761b9fc7fcSBarry Smith if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split"); 277704ba839SBarry Smith ierr = ISAllGather(ilink->is,&ilink->cis);CHKERRQ(ierr); 278704ba839SBarry Smith ilink = ilink->next; 2791b9fc7fcSBarry Smith } 2801b9fc7fcSBarry Smith } 2811b9fc7fcSBarry Smith 282704ba839SBarry Smith ilink = jac->head; 28397bbdb24SBarry Smith if (!jac->pmat) { 284cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 285cf502942SBarry Smith for (i=0; i<nsplit; i++) { 286704ba839SBarry Smith ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 287704ba839SBarry Smith ilink = ilink->next; 288cf502942SBarry Smith } 28997bbdb24SBarry Smith } else { 290cf502942SBarry Smith for (i=0; i<nsplit; i++) { 291704ba839SBarry Smith ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 292704ba839SBarry Smith ilink = ilink->next; 293cf502942SBarry Smith } 29497bbdb24SBarry Smith } 29597bbdb24SBarry Smith 29668dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 29768dd23aaSBarry Smith ierr = MatHasOperation(pc->mat,MATOP_GET_SUBMATRIX,&getsub);CHKERRQ(ierr); 2983b224e63SBarry Smith if (getsub && jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 29968dd23aaSBarry Smith ilink = jac->head; 30068dd23aaSBarry Smith if (!jac->Afield) { 30168dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 30268dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 30311755939SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,PETSC_DECIDE,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 30468dd23aaSBarry Smith ilink = ilink->next; 30568dd23aaSBarry Smith } 30668dd23aaSBarry Smith } else { 30768dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 30811755939SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,PETSC_DECIDE,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 30968dd23aaSBarry Smith ilink = ilink->next; 31068dd23aaSBarry Smith } 31168dd23aaSBarry Smith } 31268dd23aaSBarry Smith } 31368dd23aaSBarry Smith 3143b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 3153b224e63SBarry Smith IS ccis; 316e69d4d44SBarry Smith PetscInt N,nlocal,nis; 3173b224e63SBarry Smith if (nsplit != 2) SETERRQ(PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 31868dd23aaSBarry Smith 3193b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 3203b224e63SBarry Smith if (jac->schur) { 3213b224e63SBarry Smith ierr = MatGetSize(pc->mat,PETSC_NULL,&N);CHKERRQ(ierr); 3223b224e63SBarry Smith ilink = jac->head; 323edf189efSBarry Smith ierr = ISComplement(ilink->cis,0,N,&ccis);CHKERRQ(ierr); 324e69d4d44SBarry Smith ierr = ISGetLocalSize(ilink->is,&nis);CHKERRQ(ierr); 325e69d4d44SBarry Smith ierr = MatGetLocalSize(pc->mat,PETSC_NULL,&nlocal);CHKERRQ(ierr); 326e69d4d44SBarry Smith nlocal = nlocal - nis; 327e69d4d44SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,nlocal,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 3283b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 3293b224e63SBarry Smith ilink = ilink->next; 330edf189efSBarry Smith ierr = ISComplement(ilink->cis,0,N,&ccis);CHKERRQ(ierr); 331e69d4d44SBarry Smith ierr = ISGetLocalSize(ilink->is,&nis);CHKERRQ(ierr); 332e69d4d44SBarry Smith ierr = MatGetLocalSize(pc->mat,PETSC_NULL,&nlocal);CHKERRQ(ierr); 333e69d4d44SBarry Smith nlocal = nlocal - nis; 334e69d4d44SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,nlocal,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 3353b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 336*084e4875SJed Brown ierr = MatSchurComplementUpdate(jac->schur,jac->pmat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr); 337*084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 3383b224e63SBarry Smith 3393b224e63SBarry Smith } else { 3401cee3971SBarry Smith KSP ksp; 3413b224e63SBarry Smith 3423b224e63SBarry Smith /* extract the B and C matrices */ 3433b224e63SBarry Smith ierr = MatGetSize(pc->mat,PETSC_NULL,&N);CHKERRQ(ierr); 3443b224e63SBarry Smith ilink = jac->head; 345edf189efSBarry Smith ierr = ISComplement(ilink->cis,0,N,&ccis);CHKERRQ(ierr); 346e69d4d44SBarry Smith ierr = ISGetLocalSize(ilink->is,&nis);CHKERRQ(ierr); 347e69d4d44SBarry Smith ierr = MatGetLocalSize(pc->mat,PETSC_NULL,&nlocal);CHKERRQ(ierr); 348e69d4d44SBarry Smith nlocal = nlocal - nis; 349e69d4d44SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,nlocal,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 3503b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 3513b224e63SBarry Smith ilink = ilink->next; 352edf189efSBarry Smith ierr = ISComplement(ilink->cis,0,N,&ccis);CHKERRQ(ierr); 353e69d4d44SBarry Smith ierr = ISGetLocalSize(ilink->is,&nis);CHKERRQ(ierr); 354e69d4d44SBarry Smith ierr = MatGetLocalSize(pc->mat,PETSC_NULL,&nlocal);CHKERRQ(ierr); 355e69d4d44SBarry Smith nlocal = nlocal - nis; 356e69d4d44SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,nlocal,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 3573b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 358*084e4875SJed Brown /* Better would be to use 'mat[0]' (diagonal block of the real matrix) preconditioned by pmat[0] */ 359*084e4875SJed Brown ierr = MatCreateSchurComplement(jac->pmat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],&jac->schur);CHKERRQ(ierr); 3601cee3971SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 361e69d4d44SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr); 3623b224e63SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 3633b224e63SBarry Smith 3643b224e63SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr); 3651cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 366*084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 367*084e4875SJed Brown if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 368*084e4875SJed Brown PC pc; 369*084e4875SJed Brown ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 370*084e4875SJed Brown ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 371*084e4875SJed Brown /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 372e69d4d44SBarry Smith } 3733b224e63SBarry Smith ierr = KSPSetOptionsPrefix(jac->kspschur,((PetscObject)pc)->prefix);CHKERRQ(ierr); 374edf189efSBarry Smith ierr = KSPAppendOptionsPrefix(jac->kspschur,"fieldsplit_1_");CHKERRQ(ierr); 3753b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 3763b224e63SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 3773b224e63SBarry Smith 3783b224e63SBarry Smith ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr); 3793b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr); 3803b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr); 3813b224e63SBarry Smith ilink = jac->head; 3823b224e63SBarry Smith ilink->x = jac->x[0]; ilink->y = jac->y[0]; 3833b224e63SBarry Smith ilink = ilink->next; 3843b224e63SBarry Smith ilink->x = jac->x[1]; ilink->y = jac->y[1]; 3853b224e63SBarry Smith } 3863b224e63SBarry Smith } else { 38797bbdb24SBarry Smith /* set up the individual PCs */ 38897bbdb24SBarry Smith i = 0; 3895a9f2f41SSatish Balay ilink = jac->head; 3905a9f2f41SSatish Balay while (ilink) { 3915a9f2f41SSatish Balay ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr); 3923b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 3935a9f2f41SSatish Balay ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr); 3945a9f2f41SSatish Balay ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 39597bbdb24SBarry Smith i++; 3965a9f2f41SSatish Balay ilink = ilink->next; 3970971522cSBarry Smith } 39897bbdb24SBarry Smith 39997bbdb24SBarry Smith /* create work vectors for each split */ 4001b9fc7fcSBarry Smith if (!jac->x) { 40197bbdb24SBarry Smith ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 4025a9f2f41SSatish Balay ilink = jac->head; 40397bbdb24SBarry Smith for (i=0; i<nsplit; i++) { 404906ed7ccSBarry Smith Vec *vl,*vr; 405906ed7ccSBarry Smith 406906ed7ccSBarry Smith ierr = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr); 407906ed7ccSBarry Smith ilink->x = *vr; 408906ed7ccSBarry Smith ilink->y = *vl; 409906ed7ccSBarry Smith ierr = PetscFree(vr);CHKERRQ(ierr); 410906ed7ccSBarry Smith ierr = PetscFree(vl);CHKERRQ(ierr); 4115a9f2f41SSatish Balay jac->x[i] = ilink->x; 4125a9f2f41SSatish Balay jac->y[i] = ilink->y; 4135a9f2f41SSatish Balay ilink = ilink->next; 41497bbdb24SBarry Smith } 4153b224e63SBarry Smith } 4163b224e63SBarry Smith } 4173b224e63SBarry Smith 4183b224e63SBarry Smith 419d0f46423SBarry Smith if (!jac->head->sctx) { 4203b224e63SBarry Smith Vec xtmp; 4213b224e63SBarry Smith 42279416396SBarry Smith /* compute scatter contexts needed by multiplicative versions and non-default splits */ 4231b9fc7fcSBarry Smith 4245a9f2f41SSatish Balay ilink = jac->head; 4251b9fc7fcSBarry Smith ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 4261b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 427704ba839SBarry Smith ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 4285a9f2f41SSatish Balay ilink = ilink->next; 4291b9fc7fcSBarry Smith } 4301b9fc7fcSBarry Smith ierr = VecDestroy(xtmp);CHKERRQ(ierr); 4311b9fc7fcSBarry Smith } 4320971522cSBarry Smith PetscFunctionReturn(0); 4330971522cSBarry Smith } 4340971522cSBarry Smith 4355a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 436ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 437ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 4385a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 439ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 440ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 44179416396SBarry Smith 4420971522cSBarry Smith #undef __FUNCT__ 4433b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 4443b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 4453b224e63SBarry Smith { 4463b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 4473b224e63SBarry Smith PetscErrorCode ierr; 4483b224e63SBarry Smith KSP ksp; 4493b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 4503b224e63SBarry Smith 4513b224e63SBarry Smith PetscFunctionBegin; 4523b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 4533b224e63SBarry Smith 4543b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4553b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4563b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 4573b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 4583b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 4593b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4603b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4613b224e63SBarry Smith 4623b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 4633b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4643b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4653b224e63SBarry Smith 4663b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 4673b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 4683b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 4693b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4703b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4713b224e63SBarry Smith 4723b224e63SBarry Smith PetscFunctionReturn(0); 4733b224e63SBarry Smith } 4743b224e63SBarry Smith 4753b224e63SBarry Smith #undef __FUNCT__ 4760971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 4770971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 4780971522cSBarry Smith { 4790971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 4800971522cSBarry Smith PetscErrorCode ierr; 4815a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 4823e197d65SBarry Smith PetscInt bs,cnt; 4830971522cSBarry Smith 4840971522cSBarry Smith PetscFunctionBegin; 48551f519a2SBarry Smith CHKMEMQ; 486e25d487eSBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 48751f519a2SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 48851f519a2SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 48951f519a2SBarry Smith 49079416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 4911b9fc7fcSBarry Smith if (jac->defaultsplit) { 4920971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 4935a9f2f41SSatish Balay while (ilink) { 4945a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 4955a9f2f41SSatish Balay ilink = ilink->next; 4960971522cSBarry Smith } 4970971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 4981b9fc7fcSBarry Smith } else { 499efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 5005a9f2f41SSatish Balay while (ilink) { 5015a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 5025a9f2f41SSatish Balay ilink = ilink->next; 5031b9fc7fcSBarry Smith } 5041b9fc7fcSBarry Smith } 50516913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 50679416396SBarry Smith if (!jac->w1) { 50779416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 50879416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 50979416396SBarry Smith } 510efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 5115a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 5123e197d65SBarry Smith cnt = 1; 5135a9f2f41SSatish Balay while (ilink->next) { 5145a9f2f41SSatish Balay ilink = ilink->next; 5153e197d65SBarry Smith if (jac->Afield) { 5163e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 5173e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 5183e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 5193e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5203e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5213e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 5223e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5233e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5243e197d65SBarry Smith } else { 5253e197d65SBarry Smith /* compute the residual over the entire vector */ 5261ab39975SBarry Smith ierr = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr); 527efb30889SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 5285a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 52979416396SBarry Smith } 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 if (jac->Afield) { 53611755939SBarry Smith /* compute the residual only over the part of the vector needed */ 53711755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 53811755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 53911755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 54011755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 54111755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 54211755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 54311755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 54411755939SBarry Smith } else { 5451ab39975SBarry Smith ierr = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr); 54651f519a2SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 54751f519a2SBarry Smith ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 54879416396SBarry Smith } 54951f519a2SBarry Smith } 55011755939SBarry Smith } 55116913363SBarry Smith } else SETERRQ1(PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 55251f519a2SBarry Smith CHKMEMQ; 5530971522cSBarry Smith PetscFunctionReturn(0); 5540971522cSBarry Smith } 5550971522cSBarry Smith 556421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 557ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 558ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 559421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 560ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 561ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 562421e10b8SBarry Smith 563421e10b8SBarry Smith #undef __FUNCT__ 564421e10b8SBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 565421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 566421e10b8SBarry Smith { 567421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 568421e10b8SBarry Smith PetscErrorCode ierr; 569421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 570421e10b8SBarry Smith PetscInt bs; 571421e10b8SBarry Smith 572421e10b8SBarry Smith PetscFunctionBegin; 573421e10b8SBarry Smith CHKMEMQ; 574421e10b8SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 575421e10b8SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 576421e10b8SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 577421e10b8SBarry Smith 578421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 579421e10b8SBarry Smith if (jac->defaultsplit) { 580421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 581421e10b8SBarry Smith while (ilink) { 582421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 583421e10b8SBarry Smith ilink = ilink->next; 584421e10b8SBarry Smith } 585421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 586421e10b8SBarry Smith } else { 587421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 588421e10b8SBarry Smith while (ilink) { 589421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 590421e10b8SBarry Smith ilink = ilink->next; 591421e10b8SBarry Smith } 592421e10b8SBarry Smith } 593421e10b8SBarry Smith } else { 594421e10b8SBarry Smith if (!jac->w1) { 595421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 596421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 597421e10b8SBarry Smith } 598421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 599421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 600421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 601421e10b8SBarry Smith while (ilink->next) { 602421e10b8SBarry Smith ilink = ilink->next; 6039989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 604421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 605421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 606421e10b8SBarry Smith } 607421e10b8SBarry Smith while (ilink->previous) { 608421e10b8SBarry Smith ilink = ilink->previous; 6099989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 610421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 611421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 612421e10b8SBarry Smith } 613421e10b8SBarry Smith } else { 614421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 615421e10b8SBarry Smith ilink = ilink->next; 616421e10b8SBarry Smith } 617421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 618421e10b8SBarry Smith while (ilink->previous) { 619421e10b8SBarry Smith ilink = ilink->previous; 6209989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 621421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 622421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 623421e10b8SBarry Smith } 624421e10b8SBarry Smith } 625421e10b8SBarry Smith } 626421e10b8SBarry Smith CHKMEMQ; 627421e10b8SBarry Smith PetscFunctionReturn(0); 628421e10b8SBarry Smith } 629421e10b8SBarry Smith 6300971522cSBarry Smith #undef __FUNCT__ 6310971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 6320971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 6330971522cSBarry Smith { 6340971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6350971522cSBarry Smith PetscErrorCode ierr; 6365a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 6370971522cSBarry Smith 6380971522cSBarry Smith PetscFunctionBegin; 6395a9f2f41SSatish Balay while (ilink) { 6405a9f2f41SSatish Balay ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr); 6415a9f2f41SSatish Balay if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);} 6425a9f2f41SSatish Balay if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);} 6435a9f2f41SSatish Balay if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);} 644704ba839SBarry Smith if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);} 645704ba839SBarry Smith if (ilink->cis) {ierr = ISDestroy(ilink->cis);CHKERRQ(ierr);} 6465a9f2f41SSatish Balay next = ilink->next; 647704ba839SBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 648704ba839SBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 6495a9f2f41SSatish Balay ilink = next; 6500971522cSBarry Smith } 65105b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 65297bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 65368dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 65479416396SBarry Smith if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);} 65579416396SBarry Smith if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);} 6563b224e63SBarry Smith if (jac->schur) {ierr = MatDestroy(jac->schur);CHKERRQ(ierr);} 657*084e4875SJed Brown if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);} 658d0f46423SBarry Smith if (jac->kspschur) {ierr = KSPDestroy(jac->kspschur);CHKERRQ(ierr);} 6593b224e63SBarry Smith if (jac->B) {ierr = MatDestroy(jac->B);CHKERRQ(ierr);} 6603b224e63SBarry Smith if (jac->C) {ierr = MatDestroy(jac->C);CHKERRQ(ierr);} 6610971522cSBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 6620971522cSBarry Smith PetscFunctionReturn(0); 6630971522cSBarry Smith } 6640971522cSBarry Smith 6650971522cSBarry Smith #undef __FUNCT__ 6660971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 6670971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 6680971522cSBarry Smith { 6691b9fc7fcSBarry Smith PetscErrorCode ierr; 67051f519a2SBarry Smith PetscInt i = 0,nfields,*fields,bs; 671*084e4875SJed Brown PetscTruth flg; 6721b9fc7fcSBarry Smith char optionname[128]; 6739dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6743b224e63SBarry Smith PCCompositeType ctype; 6751b9fc7fcSBarry Smith 6760971522cSBarry Smith PetscFunctionBegin; 6771b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 67851f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 67951f519a2SBarry Smith if (flg) { 68051f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 68151f519a2SBarry Smith } 682704ba839SBarry Smith 6833b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 6843b224e63SBarry Smith if (flg) { 6853b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 6863b224e63SBarry Smith } 687d32f9abdSBarry Smith 688c30613efSMatthew Knepley /* Only setup fields once */ 689c30613efSMatthew Knepley if ((jac->bs > 0) && (jac->nsplits == 0)) { 690d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 691d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 69251f519a2SBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr); 6931b9fc7fcSBarry Smith while (PETSC_TRUE) { 69413f74950SBarry Smith sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 69551f519a2SBarry Smith nfields = jac->bs; 6961b9fc7fcSBarry Smith ierr = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr); 6971b9fc7fcSBarry Smith if (!flg) break; 6981b9fc7fcSBarry Smith if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 6991b9fc7fcSBarry Smith ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr); 7001b9fc7fcSBarry Smith } 70151f519a2SBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 702d32f9abdSBarry Smith } 703*084e4875SJed 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); 7041b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 7050971522cSBarry Smith PetscFunctionReturn(0); 7060971522cSBarry Smith } 7070971522cSBarry Smith 7080971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 7090971522cSBarry Smith 7100971522cSBarry Smith EXTERN_C_BEGIN 7110971522cSBarry Smith #undef __FUNCT__ 7120971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 713dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields) 7140971522cSBarry Smith { 71597bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7160971522cSBarry Smith PetscErrorCode ierr; 7175a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 71869a612a9SBarry Smith char prefix[128]; 71951f519a2SBarry Smith PetscInt i; 7200971522cSBarry Smith 7210971522cSBarry Smith PetscFunctionBegin; 7220971522cSBarry Smith if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested"); 72351f519a2SBarry Smith for (i=0; i<n; i++) { 72451f519a2SBarry Smith if (fields[i] >= jac->bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs); 72551f519a2SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 72651f519a2SBarry Smith } 727704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 728704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 7295a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 7305a9f2f41SSatish Balay ilink->nfields = n; 7315a9f2f41SSatish Balay ilink->next = PETSC_NULL; 7327adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 7331cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 7345a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 73569a612a9SBarry Smith 7367adad957SLisandro Dalcin if (((PetscObject)pc)->prefix) { 7377adad957SLisandro Dalcin sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 73869a612a9SBarry Smith } else { 73913f74950SBarry Smith sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 74069a612a9SBarry Smith } 7415a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 7420971522cSBarry Smith 7430971522cSBarry Smith if (!next) { 7445a9f2f41SSatish Balay jac->head = ilink; 74551f519a2SBarry Smith ilink->previous = PETSC_NULL; 7460971522cSBarry Smith } else { 7470971522cSBarry Smith while (next->next) { 7480971522cSBarry Smith next = next->next; 7490971522cSBarry Smith } 7505a9f2f41SSatish Balay next->next = ilink; 75151f519a2SBarry Smith ilink->previous = next; 7520971522cSBarry Smith } 7530971522cSBarry Smith jac->nsplits++; 7540971522cSBarry Smith PetscFunctionReturn(0); 7550971522cSBarry Smith } 7560971522cSBarry Smith EXTERN_C_END 7570971522cSBarry Smith 758e69d4d44SBarry Smith EXTERN_C_BEGIN 759e69d4d44SBarry Smith #undef __FUNCT__ 760e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 761e69d4d44SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 762e69d4d44SBarry Smith { 763e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 764e69d4d44SBarry Smith PetscErrorCode ierr; 765e69d4d44SBarry Smith 766e69d4d44SBarry Smith PetscFunctionBegin; 767e69d4d44SBarry Smith ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 768e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 769e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 770*084e4875SJed Brown *n = jac->nsplits; 771e69d4d44SBarry Smith PetscFunctionReturn(0); 772e69d4d44SBarry Smith } 773e69d4d44SBarry Smith EXTERN_C_END 7740971522cSBarry Smith 7750971522cSBarry Smith EXTERN_C_BEGIN 7760971522cSBarry Smith #undef __FUNCT__ 77769a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 778dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 7790971522cSBarry Smith { 7800971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7810971522cSBarry Smith PetscErrorCode ierr; 7820971522cSBarry Smith PetscInt cnt = 0; 7835a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 7840971522cSBarry Smith 7850971522cSBarry Smith PetscFunctionBegin; 78669a612a9SBarry Smith ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 7875a9f2f41SSatish Balay while (ilink) { 7885a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 7895a9f2f41SSatish Balay ilink = ilink->next; 7900971522cSBarry Smith } 7910971522cSBarry Smith if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits); 7920971522cSBarry Smith *n = jac->nsplits; 7930971522cSBarry Smith PetscFunctionReturn(0); 7940971522cSBarry Smith } 7950971522cSBarry Smith EXTERN_C_END 7960971522cSBarry Smith 797704ba839SBarry Smith EXTERN_C_BEGIN 798704ba839SBarry Smith #undef __FUNCT__ 799704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 800704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,IS is) 801704ba839SBarry Smith { 802704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 803704ba839SBarry Smith PetscErrorCode ierr; 804704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 805704ba839SBarry Smith char prefix[128]; 806704ba839SBarry Smith 807704ba839SBarry Smith PetscFunctionBegin; 80816913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 8091ab39975SBarry Smith ilink->is = is; 8101ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 811704ba839SBarry Smith ilink->next = PETSC_NULL; 812704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 8131cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 814704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 815704ba839SBarry Smith 816704ba839SBarry Smith if (((PetscObject)pc)->prefix) { 817704ba839SBarry Smith sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 818704ba839SBarry Smith } else { 819704ba839SBarry Smith sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 820704ba839SBarry Smith } 821704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 822704ba839SBarry Smith 823704ba839SBarry Smith if (!next) { 824704ba839SBarry Smith jac->head = ilink; 825704ba839SBarry Smith ilink->previous = PETSC_NULL; 826704ba839SBarry Smith } else { 827704ba839SBarry Smith while (next->next) { 828704ba839SBarry Smith next = next->next; 829704ba839SBarry Smith } 830704ba839SBarry Smith next->next = ilink; 831704ba839SBarry Smith ilink->previous = next; 832704ba839SBarry Smith } 833704ba839SBarry Smith jac->nsplits++; 834704ba839SBarry Smith 835704ba839SBarry Smith PetscFunctionReturn(0); 836704ba839SBarry Smith } 837704ba839SBarry Smith EXTERN_C_END 838704ba839SBarry Smith 8390971522cSBarry Smith #undef __FUNCT__ 8400971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 8410971522cSBarry Smith /*@ 8420971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 8430971522cSBarry Smith 8440971522cSBarry Smith Collective on PC 8450971522cSBarry Smith 8460971522cSBarry Smith Input Parameters: 8470971522cSBarry Smith + pc - the preconditioner context 8480971522cSBarry Smith . n - the number of fields in this split 8490971522cSBarry Smith . fields - the fields in this split 8500971522cSBarry Smith 8510971522cSBarry Smith Level: intermediate 8520971522cSBarry Smith 853d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 854d32f9abdSBarry Smith 855d32f9abdSBarry Smith The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block 856d32f9abdSBarry 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 857d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 858d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 859d32f9abdSBarry Smith 860d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 8610971522cSBarry Smith 8620971522cSBarry Smith @*/ 863dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields) 8640971522cSBarry Smith { 8650971522cSBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *); 8660971522cSBarry Smith 8670971522cSBarry Smith PetscFunctionBegin; 8680971522cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 8690971522cSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr); 8700971522cSBarry Smith if (f) { 8710971522cSBarry Smith ierr = (*f)(pc,n,fields);CHKERRQ(ierr); 8720971522cSBarry Smith } 8730971522cSBarry Smith PetscFunctionReturn(0); 8740971522cSBarry Smith } 8750971522cSBarry Smith 8760971522cSBarry Smith #undef __FUNCT__ 877704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 878704ba839SBarry Smith /*@ 879704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 880704ba839SBarry Smith 881704ba839SBarry Smith Collective on PC 882704ba839SBarry Smith 883704ba839SBarry Smith Input Parameters: 884704ba839SBarry Smith + pc - the preconditioner context 885704ba839SBarry Smith . is - the index set that defines the vector elements in this field 886704ba839SBarry Smith 887d32f9abdSBarry Smith 888d32f9abdSBarry Smith Notes: Use PCFieldSplitSetFields(), for fields defined by strided types. 889d32f9abdSBarry Smith 890704ba839SBarry Smith Level: intermediate 891704ba839SBarry Smith 892704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 893704ba839SBarry Smith 894704ba839SBarry Smith @*/ 895704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,IS is) 896704ba839SBarry Smith { 897704ba839SBarry Smith PetscErrorCode ierr,(*f)(PC,IS); 898704ba839SBarry Smith 899704ba839SBarry Smith PetscFunctionBegin; 900704ba839SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 901704ba839SBarry Smith PetscValidHeaderSpecific(is,IS_COOKIE,1); 902704ba839SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr); 903704ba839SBarry Smith if (f) { 904704ba839SBarry Smith ierr = (*f)(pc,is);CHKERRQ(ierr); 905704ba839SBarry Smith } 906704ba839SBarry Smith PetscFunctionReturn(0); 907704ba839SBarry Smith } 908704ba839SBarry Smith 909704ba839SBarry Smith #undef __FUNCT__ 91051f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 91151f519a2SBarry Smith /*@ 91251f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 91351f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 91451f519a2SBarry Smith 91551f519a2SBarry Smith Collective on PC 91651f519a2SBarry Smith 91751f519a2SBarry Smith Input Parameters: 91851f519a2SBarry Smith + pc - the preconditioner context 91951f519a2SBarry Smith - bs - the block size 92051f519a2SBarry Smith 92151f519a2SBarry Smith Level: intermediate 92251f519a2SBarry Smith 92351f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 92451f519a2SBarry Smith 92551f519a2SBarry Smith @*/ 92651f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 92751f519a2SBarry Smith { 92851f519a2SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt); 92951f519a2SBarry Smith 93051f519a2SBarry Smith PetscFunctionBegin; 93151f519a2SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 93251f519a2SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr); 93351f519a2SBarry Smith if (f) { 93451f519a2SBarry Smith ierr = (*f)(pc,bs);CHKERRQ(ierr); 93551f519a2SBarry Smith } 93651f519a2SBarry Smith PetscFunctionReturn(0); 93751f519a2SBarry Smith } 93851f519a2SBarry Smith 93951f519a2SBarry Smith #undef __FUNCT__ 94069a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 9410971522cSBarry Smith /*@C 94269a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 9430971522cSBarry Smith 94469a612a9SBarry Smith Collective on KSP 9450971522cSBarry Smith 9460971522cSBarry Smith Input Parameter: 9470971522cSBarry Smith . pc - the preconditioner context 9480971522cSBarry Smith 9490971522cSBarry Smith Output Parameters: 9500971522cSBarry Smith + n - the number of split 95169a612a9SBarry Smith - pc - the array of KSP contexts 9520971522cSBarry Smith 9530971522cSBarry Smith Note: 954d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 955d32f9abdSBarry Smith (not the KSP just the array that contains them). 9560971522cSBarry Smith 95769a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 9580971522cSBarry Smith 9590971522cSBarry Smith Level: advanced 9600971522cSBarry Smith 9610971522cSBarry Smith .seealso: PCFIELDSPLIT 9620971522cSBarry Smith @*/ 963dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 9640971522cSBarry Smith { 96569a612a9SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **); 9660971522cSBarry Smith 9670971522cSBarry Smith PetscFunctionBegin; 9680971522cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 9690971522cSBarry Smith PetscValidIntPointer(n,2); 97069a612a9SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 9710971522cSBarry Smith if (f) { 97269a612a9SBarry Smith ierr = (*f)(pc,n,subksp);CHKERRQ(ierr); 9730971522cSBarry Smith } else { 97469a612a9SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC"); 9750971522cSBarry Smith } 9760971522cSBarry Smith PetscFunctionReturn(0); 9770971522cSBarry Smith } 9780971522cSBarry Smith 979e69d4d44SBarry Smith #undef __FUNCT__ 980e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 981e69d4d44SBarry Smith /*@ 982e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 983e69d4d44SBarry Smith D matrix. Otherwise no preconditioner is used. 984e69d4d44SBarry Smith 985e69d4d44SBarry Smith Collective on PC 986e69d4d44SBarry Smith 987e69d4d44SBarry Smith Input Parameters: 988e69d4d44SBarry Smith + pc - the preconditioner context 989*084e4875SJed Brown . ptype - which matrix to use for preconditioning the Schur complement 990*084e4875SJed Brown - userpre - matrix to use for preconditioning, or PETSC_NULL 991*084e4875SJed Brown 992*084e4875SJed Brown Notes: 993*084e4875SJed Brown The default is to use the block on the diagonal of the preconditioning matrix. This is D, in the (1,1) position. 994*084e4875SJed Brown There are currently no preconditioners that work directly with the Schur complement so setting 995*084e4875SJed Brown PC_FIELDSPLIT_SCHUR_PRE_SELF is observationally equivalent to -fieldsplit_1_pc_type none. 996e69d4d44SBarry Smith 997e69d4d44SBarry Smith Options Database: 998*084e4875SJed Brown . -pc_fieldsplit_schur_precondition <self,user,diag> default is diag 999e69d4d44SBarry Smith 1000e69d4d44SBarry Smith Level: intermediate 1001e69d4d44SBarry Smith 1002*084e4875SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType 1003e69d4d44SBarry Smith 1004e69d4d44SBarry Smith @*/ 1005*084e4875SJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1006e69d4d44SBarry Smith { 1007*084e4875SJed Brown PetscErrorCode ierr,(*f)(PC,PCFieldSplitSchurPreType,Mat); 1008e69d4d44SBarry Smith 1009e69d4d44SBarry Smith PetscFunctionBegin; 1010e69d4d44SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1011e69d4d44SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",(void (**)(void))&f);CHKERRQ(ierr); 1012e69d4d44SBarry Smith if (f) { 1013*084e4875SJed Brown ierr = (*f)(pc,ptype,pre);CHKERRQ(ierr); 1014e69d4d44SBarry Smith } 1015e69d4d44SBarry Smith PetscFunctionReturn(0); 1016e69d4d44SBarry Smith } 1017e69d4d44SBarry Smith 1018e69d4d44SBarry Smith EXTERN_C_BEGIN 1019e69d4d44SBarry Smith #undef __FUNCT__ 1020e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 1021*084e4875SJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1022e69d4d44SBarry Smith { 1023e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1024*084e4875SJed Brown PetscErrorCode ierr; 1025e69d4d44SBarry Smith 1026e69d4d44SBarry Smith PetscFunctionBegin; 1027*084e4875SJed Brown jac->schurpre = ptype; 1028*084e4875SJed Brown if (pre) { 1029*084e4875SJed Brown if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);} 1030*084e4875SJed Brown jac->schur_user = pre; 1031*084e4875SJed Brown ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1032*084e4875SJed Brown } 1033e69d4d44SBarry Smith PetscFunctionReturn(0); 1034e69d4d44SBarry Smith } 1035e69d4d44SBarry Smith EXTERN_C_END 1036e69d4d44SBarry Smith 103730ad9308SMatthew Knepley #undef __FUNCT__ 103830ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 103930ad9308SMatthew Knepley /*@C 104030ad9308SMatthew Knepley PCFieldSplitGetSchurBlocks - Gets the all matrix blocks for the Schur complement 104130ad9308SMatthew Knepley 104230ad9308SMatthew Knepley Collective on KSP 104330ad9308SMatthew Knepley 104430ad9308SMatthew Knepley Input Parameter: 104530ad9308SMatthew Knepley . pc - the preconditioner context 104630ad9308SMatthew Knepley 104730ad9308SMatthew Knepley Output Parameters: 104830ad9308SMatthew Knepley + A - the (0,0) block 104930ad9308SMatthew Knepley . B - the (0,1) block 105030ad9308SMatthew Knepley . C - the (1,0) block 105130ad9308SMatthew Knepley - D - the (1,1) block 105230ad9308SMatthew Knepley 105330ad9308SMatthew Knepley Level: advanced 105430ad9308SMatthew Knepley 105530ad9308SMatthew Knepley .seealso: PCFIELDSPLIT 105630ad9308SMatthew Knepley @*/ 105730ad9308SMatthew Knepley PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSchurBlocks(PC pc,Mat *A,Mat *B,Mat *C, Mat *D) 105830ad9308SMatthew Knepley { 105930ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 106030ad9308SMatthew Knepley 106130ad9308SMatthew Knepley PetscFunctionBegin; 106230ad9308SMatthew Knepley PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1063cf8ad460SMatthew Knepley if (jac->type != PC_COMPOSITE_SCHUR) {SETERRQ(PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");} 106430ad9308SMatthew Knepley if (A) *A = jac->pmat[0]; 106530ad9308SMatthew Knepley if (B) *B = jac->B; 106630ad9308SMatthew Knepley if (C) *C = jac->C; 106730ad9308SMatthew Knepley if (D) *D = jac->pmat[1]; 106830ad9308SMatthew Knepley PetscFunctionReturn(0); 106930ad9308SMatthew Knepley } 107030ad9308SMatthew Knepley 107179416396SBarry Smith EXTERN_C_BEGIN 107279416396SBarry Smith #undef __FUNCT__ 107379416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 1074dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 107579416396SBarry Smith { 107679416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1077e69d4d44SBarry Smith PetscErrorCode ierr; 107879416396SBarry Smith 107979416396SBarry Smith PetscFunctionBegin; 108079416396SBarry Smith jac->type = type; 10813b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 10823b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 10833b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 1084e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1085e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1086e69d4d44SBarry Smith 10873b224e63SBarry Smith } else { 10883b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 10893b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 1090e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 10919e7d6b0aSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 10923b224e63SBarry Smith } 109379416396SBarry Smith PetscFunctionReturn(0); 109479416396SBarry Smith } 109579416396SBarry Smith EXTERN_C_END 109679416396SBarry Smith 109751f519a2SBarry Smith EXTERN_C_BEGIN 109851f519a2SBarry Smith #undef __FUNCT__ 109951f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 110051f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 110151f519a2SBarry Smith { 110251f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 110351f519a2SBarry Smith 110451f519a2SBarry Smith PetscFunctionBegin; 110551f519a2SBarry Smith if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 110651f519a2SBarry 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); 110751f519a2SBarry Smith jac->bs = bs; 110851f519a2SBarry Smith PetscFunctionReturn(0); 110951f519a2SBarry Smith } 111051f519a2SBarry Smith EXTERN_C_END 111151f519a2SBarry Smith 111279416396SBarry Smith #undef __FUNCT__ 111379416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1114bc08b0f1SBarry Smith /*@ 111579416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 111679416396SBarry Smith 111779416396SBarry Smith Collective on PC 111879416396SBarry Smith 111979416396SBarry Smith Input Parameter: 112079416396SBarry Smith . pc - the preconditioner context 112181540f2fSBarry Smith . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 112279416396SBarry Smith 112379416396SBarry Smith Options Database Key: 1124a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 112579416396SBarry Smith 112679416396SBarry Smith Level: Developer 112779416396SBarry Smith 112879416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 112979416396SBarry Smith 113079416396SBarry Smith .seealso: PCCompositeSetType() 113179416396SBarry Smith 113279416396SBarry Smith @*/ 1133dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type) 113479416396SBarry Smith { 113579416396SBarry Smith PetscErrorCode ierr,(*f)(PC,PCCompositeType); 113679416396SBarry Smith 113779416396SBarry Smith PetscFunctionBegin; 113879416396SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 113979416396SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 114079416396SBarry Smith if (f) { 114179416396SBarry Smith ierr = (*f)(pc,type);CHKERRQ(ierr); 114279416396SBarry Smith } 114379416396SBarry Smith PetscFunctionReturn(0); 114479416396SBarry Smith } 114579416396SBarry Smith 11460971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 11470971522cSBarry Smith /*MC 1148a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 11490971522cSBarry Smith fields or groups of fields 11500971522cSBarry Smith 11510971522cSBarry Smith 1152edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1153edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 11540971522cSBarry Smith 1155a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 115669a612a9SBarry Smith and set the options directly on the resulting KSP object 11570971522cSBarry Smith 11580971522cSBarry Smith Level: intermediate 11590971522cSBarry Smith 116079416396SBarry Smith Options Database Keys: 116181540f2fSBarry Smith + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 116281540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 116381540f2fSBarry Smith been supplied explicitly by -pc_fieldsplit_%d_fields 116481540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 116581540f2fSBarry Smith . -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative> 1166e69d4d44SBarry Smith . -pc_fieldsplit_schur_precondition <true,false> default is true 116779416396SBarry Smith 1168edf189efSBarry Smith - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 11693b224e63SBarry Smith for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 11703b224e63SBarry Smith 11713b224e63SBarry Smith 1172d32f9abdSBarry Smith Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1173d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1174d32f9abdSBarry Smith 1175d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1176d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1177d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1178d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1179d32f9abdSBarry Smith 1180d32f9abdSBarry Smith Currently for the multiplicative version, the updated residual needed for the next field 1181d32f9abdSBarry Smith solve is computed via a matrix vector product over the entire array. An optimization would be 1182d32f9abdSBarry Smith to update the residual only for the part of the right hand side associated with the next field 1183d32f9abdSBarry Smith solve. (This would involve more MatGetSubMatrix() calls or some other mechanism to compute the 1184d32f9abdSBarry Smith part of the matrix needed to just update part of the residual). 1185d32f9abdSBarry Smith 1186e69d4d44SBarry Smith For the Schur complement preconditioner if J = ( A B ) 1187e69d4d44SBarry Smith ( C D ) 1188e69d4d44SBarry Smith the preconditioner is 1189e69d4d44SBarry Smith (I -B inv(A)) ( inv(A) 0 ) (I 0 ) 1190e69d4d44SBarry Smith (0 I ) ( 0 inv(S) ) (-C inv(A) I ) 1191edf189efSBarry Smith where the action of inv(A) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1192e69d4d44SBarry Smith inv(S) is computed using the KSP solver with prefix -schur_. For PCFieldSplitGetKSP() when field number is 1193edf189efSBarry Smith 0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1194e69d4d44SBarry Smith D is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1195e69d4d44SBarry Smith option. 1196e69d4d44SBarry Smith 1197edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1198edf189efSBarry Smith is used automatically for a second block. 1199edf189efSBarry Smith 1200a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 12010971522cSBarry Smith 1202a541d17aSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners 1203e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 12040971522cSBarry Smith M*/ 12050971522cSBarry Smith 12060971522cSBarry Smith EXTERN_C_BEGIN 12070971522cSBarry Smith #undef __FUNCT__ 12080971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 1209dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc) 12100971522cSBarry Smith { 12110971522cSBarry Smith PetscErrorCode ierr; 12120971522cSBarry Smith PC_FieldSplit *jac; 12130971522cSBarry Smith 12140971522cSBarry Smith PetscFunctionBegin; 121538f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 12160971522cSBarry Smith jac->bs = -1; 12170971522cSBarry Smith jac->nsplits = 0; 12183e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1219*084e4875SJed Brown jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_DIAG; 122051f519a2SBarry Smith 12210971522cSBarry Smith pc->data = (void*)jac; 12220971522cSBarry Smith 12230971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1224421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 12250971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 12260971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 12270971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 12280971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 12290971522cSBarry Smith pc->ops->applyrichardson = 0; 12300971522cSBarry Smith 123169a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 123269a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 12330971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 12340971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1235704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1236704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 123779416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 123879416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 123951f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 124051f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 12410971522cSBarry Smith PetscFunctionReturn(0); 12420971522cSBarry Smith } 12430971522cSBarry Smith EXTERN_C_END 12440971522cSBarry Smith 12450971522cSBarry Smith 1246a541d17aSBarry Smith /*MC 1247a541d17aSBarry Smith Block_Preconditioners - PETSc provides support for a variety of "block preconditioners", this provides an 1248a541d17aSBarry Smith overview of these methods. 1249a541d17aSBarry Smith 1250a541d17aSBarry Smith Consider the solution to ( A B ) (x_1) = (b_1) 1251a541d17aSBarry Smith ( C D ) (x_2) (b_2) 1252a541d17aSBarry Smith 1253a541d17aSBarry Smith Important special cases, the Stokes equation: C = B' and D = 0 (A B) (x_1) = (b_1) 1254a541d17aSBarry Smith B' 0) (x_2) (b_2) 1255a541d17aSBarry Smith 1256a541d17aSBarry Smith One of the goals of the PCFieldSplit preconditioner in PETSc is to provide a variety of preconditioners 1257a541d17aSBarry Smith for this block system. 1258a541d17aSBarry Smith 1259a541d17aSBarry Smith Consider an additional matrix (Ap Bp) 1260a541d17aSBarry Smith (Cp Dp) where some or all of the entries may be the same as 1261a541d17aSBarry Smith in the original matrix (for example Ap == A). 1262a541d17aSBarry Smith 1263a541d17aSBarry Smith In the following, A^ denotes the approximate application of the inverse of A, possibly using Ap in the 1264a541d17aSBarry Smith approximation. In PETSc this simply means one has called KSPSetOperators(ksp,A,Ap,...) or KSPSetOperators(ksp,Ap,Ap,...) 1265a541d17aSBarry Smith 1266a541d17aSBarry Smith Block Jacobi: x_1 = A^ b_1 1267a541d17aSBarry Smith x_2 = D^ b_2 1268a541d17aSBarry Smith 1269a541d17aSBarry Smith Lower block Gauss-Seidel: x_1 = A^ b_1 1270a541d17aSBarry Smith x_2 = D^ (b_2 - C x_1) variant x_2 = D^ (b_2 - Cp x_1) 1271a541d17aSBarry Smith 1272a541d17aSBarry 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) 1273a541d17aSBarry Smith Interestingly this form is not actually a symmetric matrix, the symmetric version is 1274a541d17aSBarry Smith x_1 = A^(b_1 - B x_2) variant x_1 = A^(b_1 - Bp x_2) 1275a541d17aSBarry Smith 12760bc0a719SBarry Smith Level: intermediate 12770bc0a719SBarry Smith 1278a541d17aSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCFIELDSPLIT 1279a541d17aSBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 1280a541d17aSBarry Smith M*/ 1281