1dba47a55SKris Buschelman #define PETSCKSP_DLL 2dba47a55SKris Buschelman 36356e834SBarry Smith #include "private/pcimpl.h" /*I "petscpc.h" I*/ 40971522cSBarry Smith 5084e4875SJed Brown const char *PCFieldSplitSchurPreTypes[] = {"self","diag","user","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0}; 6084e4875SJed Brown 70971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink; 80971522cSBarry Smith struct _PC_FieldSplitLink { 969a612a9SBarry Smith KSP ksp; 100971522cSBarry Smith Vec x,y; 11db4c96c1SJed Brown char *splitname; 120971522cSBarry Smith PetscInt nfields; 130971522cSBarry Smith PetscInt *fields; 141b9fc7fcSBarry Smith VecScatter sctx; 154aa3045dSJed Brown IS is; 1651f519a2SBarry Smith PC_FieldSplitLink next,previous; 170971522cSBarry Smith }; 180971522cSBarry Smith 190971522cSBarry Smith typedef struct { 2068dd23aaSBarry Smith PCCompositeType type; 21a4efd8eaSMatthew Knepley PetscTruth defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */ 226c924f48SJed Brown PetscTruth splitdefined; /* Flag is set after the splits have been defined, to prevent more splits from being added */ 23519d70e2SJed Brown PetscTruth realdiagonal; /* Flag to use the diagonal blocks of mat preconditioned by pmat, instead of just pmat */ 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; 27519d70e2SJed Brown Mat *mat; /* The diagonal block for each split */ 28519d70e2SJed Brown Mat *pmat; /* The preconditioning 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 */ 35084e4875SJed Brown Mat schur_user; /* User-provided preconditioning matrix for the Schur complement */ 36084e4875SJed 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 */ 46084e4875SJed Brown 47e6cab6aaSJed Brown /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it. This way the 48084e4875SJed Brown * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */ 49084e4875SJed Brown static Mat FieldSplitSchurPre(PC_FieldSplit *jac) 50084e4875SJed Brown { 51084e4875SJed Brown switch (jac->schurpre) { 52084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur; 53084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1]; 54084e4875SJed Brown case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */ 55084e4875SJed Brown default: 56084e4875SJed Brown return jac->schur_user ? jac->schur_user : jac->pmat[1]; 57084e4875SJed Brown } 58084e4875SJed Brown } 59084e4875SJed Brown 60084e4875SJed 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; 722692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&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 { 9765e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,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; 1142692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&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 = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 13812cae6f2SJed Brown if (jac->schur) { 13912cae6f2SJed Brown ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 1403b224e63SBarry Smith ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 14112cae6f2SJed Brown } else { 14212cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 14312cae6f2SJed Brown } 1443b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1453b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = D - C inv(A) B \n");CHKERRQ(ierr); 1463b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 14712cae6f2SJed Brown if (jac->kspschur) { 1483b224e63SBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 14912cae6f2SJed Brown } else { 15012cae6f2SJed Brown ierr = PetscViewerASCIIPrintf(viewer," not yet available\n");CHKERRQ(ierr); 15112cae6f2SJed Brown } 1523b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1533b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1543b224e63SBarry Smith } else { 15565e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1563b224e63SBarry Smith } 1573b224e63SBarry Smith PetscFunctionReturn(0); 1583b224e63SBarry Smith } 1593b224e63SBarry Smith 1603b224e63SBarry Smith #undef __FUNCT__ 1616c924f48SJed Brown #define __FUNCT__ "PCFieldSplitSetRuntimeSplits_Private" 1626c924f48SJed Brown /* Precondition: jac->bs is set to a meaningful value */ 1636c924f48SJed Brown static PetscErrorCode PCFieldSplitSetRuntimeSplits_Private(PC pc) 1646c924f48SJed Brown { 1656c924f48SJed Brown PetscErrorCode ierr; 1666c924f48SJed Brown PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1676c924f48SJed Brown PetscInt i,nfields,*ifields; 1686c924f48SJed Brown PetscTruth flg; 1696c924f48SJed Brown char optionname[128],splitname[8]; 1706c924f48SJed Brown 1716c924f48SJed Brown PetscFunctionBegin; 1726c924f48SJed Brown ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 1736c924f48SJed Brown for (i=0,flg=PETSC_TRUE; ; i++) { 1746c924f48SJed Brown ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 1756c924f48SJed Brown ierr = PetscSNPrintf(optionname,sizeof optionname,"-pc_fieldsplit_%D_fields",i);CHKERRQ(ierr); 1766c924f48SJed Brown nfields = jac->bs; 1776c924f48SJed Brown ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg);CHKERRQ(ierr); 1786c924f48SJed Brown if (!flg) break; 1796c924f48SJed Brown if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields"); 1806c924f48SJed Brown ierr = PCFieldSplitSetFields(pc,splitname,nfields,ifields);CHKERRQ(ierr); 1816c924f48SJed Brown } 1826c924f48SJed Brown if (i > 0) { 1836c924f48SJed Brown /* Makes command-line setting of splits take precedence over setting them in code. 1846c924f48SJed Brown Otherwise subsequent calls to PCFieldSplitSetIS() or PCFieldSplitSetFields() would 1856c924f48SJed Brown create new splits, which would probably not be what the user wanted. */ 1866c924f48SJed Brown jac->splitdefined = PETSC_TRUE; 1876c924f48SJed Brown } 1886c924f48SJed Brown ierr = PetscFree(ifields);CHKERRQ(ierr); 1896c924f48SJed Brown PetscFunctionReturn(0); 1906c924f48SJed Brown } 1916c924f48SJed Brown 1926c924f48SJed Brown #undef __FUNCT__ 19369a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 19469a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 1950971522cSBarry Smith { 1960971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1970971522cSBarry Smith PetscErrorCode ierr; 1985a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 1996c924f48SJed Brown PetscTruth flg = PETSC_FALSE; 2006c924f48SJed Brown PetscInt i; 2010971522cSBarry Smith 2020971522cSBarry Smith PetscFunctionBegin; 203d32f9abdSBarry Smith if (!ilink) { 204d32f9abdSBarry Smith 205521d7252SBarry Smith if (jac->bs <= 0) { 206704ba839SBarry Smith if (pc->pmat) { 207521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 208704ba839SBarry Smith } else { 209704ba839SBarry Smith jac->bs = 1; 210704ba839SBarry Smith } 211521d7252SBarry Smith } 212d32f9abdSBarry Smith 213d32f9abdSBarry Smith ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 214d32f9abdSBarry Smith if (!flg) { 215d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 216d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 2176c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 2186c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 219d32f9abdSBarry Smith } 2206c924f48SJed Brown if (flg || !jac->splitdefined) { 221d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 222db4c96c1SJed Brown for (i=0; i<jac->bs; i++) { 2236c924f48SJed Brown char splitname[8]; 2246c924f48SJed Brown ierr = PetscSNPrintf(splitname,sizeof splitname,"%D",i);CHKERRQ(ierr); 225db4c96c1SJed Brown ierr = PCFieldSplitSetFields(pc,splitname,1,&i);CHKERRQ(ierr); 22679416396SBarry Smith } 22797bbdb24SBarry Smith jac->defaultsplit = PETSC_TRUE; 228521d7252SBarry Smith } 229edf189efSBarry Smith } else if (jac->nsplits == 1) { 230edf189efSBarry Smith if (ilink->is) { 231edf189efSBarry Smith IS is2; 232edf189efSBarry Smith PetscInt nmin,nmax; 233edf189efSBarry Smith 234edf189efSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 235edf189efSBarry Smith ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 236db4c96c1SJed Brown ierr = PCFieldSplitSetIS(pc,"1",is2);CHKERRQ(ierr); 237edf189efSBarry Smith ierr = ISDestroy(is2);CHKERRQ(ierr); 238e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 239edf189efSBarry Smith } 240e7e72b3dSBarry Smith if (jac->nsplits < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields"); 24169a612a9SBarry Smith PetscFunctionReturn(0); 24269a612a9SBarry Smith } 24369a612a9SBarry Smith 24469a612a9SBarry Smith 24569a612a9SBarry Smith #undef __FUNCT__ 24669a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 24769a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 24869a612a9SBarry Smith { 24969a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 25069a612a9SBarry Smith PetscErrorCode ierr; 2515a9f2f41SSatish Balay PC_FieldSplitLink ilink; 252cf502942SBarry Smith PetscInt i,nsplit,ccsize; 25369a612a9SBarry Smith MatStructure flag = pc->flag; 2546c8605c2SJed Brown PetscTruth sorted; 25569a612a9SBarry Smith 25669a612a9SBarry Smith PetscFunctionBegin; 25769a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 25897bbdb24SBarry Smith nsplit = jac->nsplits; 2595a9f2f41SSatish Balay ilink = jac->head; 26097bbdb24SBarry Smith 26197bbdb24SBarry Smith /* get the matrices for each split */ 262704ba839SBarry Smith if (!jac->issetup) { 2631b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 26497bbdb24SBarry Smith 265704ba839SBarry Smith jac->issetup = PETSC_TRUE; 266704ba839SBarry Smith 267704ba839SBarry Smith /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 26851f519a2SBarry Smith bs = jac->bs; 26997bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 270cf502942SBarry Smith ierr = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr); 2711b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 2721b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 2731b9fc7fcSBarry Smith if (jac->defaultsplit) { 274704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 275704ba839SBarry Smith } else if (!ilink->is) { 276ccb205f8SBarry Smith if (ilink->nfields > 1) { 2775a9f2f41SSatish Balay PetscInt *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields; 2785a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 2791b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 2801b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 2811b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 28297bbdb24SBarry Smith } 28397bbdb24SBarry Smith } 284704ba839SBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr); 285ccb205f8SBarry Smith ierr = PetscFree(ii);CHKERRQ(ierr); 286ccb205f8SBarry Smith } else { 287704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 288ccb205f8SBarry Smith } 2893e197d65SBarry Smith } 290704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 291e32f2f54SBarry Smith if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split"); 292704ba839SBarry Smith ilink = ilink->next; 2931b9fc7fcSBarry Smith } 2941b9fc7fcSBarry Smith } 2951b9fc7fcSBarry Smith 296704ba839SBarry Smith ilink = jac->head; 29797bbdb24SBarry Smith if (!jac->pmat) { 298cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 299cf502942SBarry Smith for (i=0; i<nsplit; i++) { 3004aa3045dSJed Brown ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 301704ba839SBarry Smith ilink = ilink->next; 302cf502942SBarry Smith } 30397bbdb24SBarry Smith } else { 304cf502942SBarry Smith for (i=0; i<nsplit; i++) { 3054aa3045dSJed Brown ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 306704ba839SBarry Smith ilink = ilink->next; 307cf502942SBarry Smith } 30897bbdb24SBarry Smith } 309519d70e2SJed Brown if (jac->realdiagonal) { 310519d70e2SJed Brown ilink = jac->head; 311519d70e2SJed Brown if (!jac->mat) { 312519d70e2SJed Brown ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr); 313519d70e2SJed Brown for (i=0; i<nsplit; i++) { 314519d70e2SJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 315519d70e2SJed Brown ilink = ilink->next; 316519d70e2SJed Brown } 317519d70e2SJed Brown } else { 318519d70e2SJed Brown for (i=0; i<nsplit; i++) { 319519d70e2SJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr); 320519d70e2SJed Brown ilink = ilink->next; 321519d70e2SJed Brown } 322519d70e2SJed Brown } 323519d70e2SJed Brown } else { 324519d70e2SJed Brown jac->mat = jac->pmat; 325519d70e2SJed Brown } 32697bbdb24SBarry Smith 3276c8605c2SJed Brown if (jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 32868dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 32968dd23aaSBarry Smith ilink = jac->head; 33068dd23aaSBarry Smith if (!jac->Afield) { 33168dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 33268dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 3334aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 33468dd23aaSBarry Smith ilink = ilink->next; 33568dd23aaSBarry Smith } 33668dd23aaSBarry Smith } else { 33768dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 3384aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 33968dd23aaSBarry Smith ilink = ilink->next; 34068dd23aaSBarry Smith } 34168dd23aaSBarry Smith } 34268dd23aaSBarry Smith } 34368dd23aaSBarry Smith 3443b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 3453b224e63SBarry Smith IS ccis; 3464aa3045dSJed Brown PetscInt rstart,rend; 347e7e72b3dSBarry Smith if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 34868dd23aaSBarry Smith 349e6cab6aaSJed Brown /* When extracting off-diagonal submatrices, we take complements from this range */ 350e6cab6aaSJed Brown ierr = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr); 351e6cab6aaSJed Brown 3523b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 3533b224e63SBarry Smith if (jac->schur) { 3543b224e63SBarry Smith ilink = jac->head; 3554aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 3564aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 3573b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 3583b224e63SBarry Smith ilink = ilink->next; 3594aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 3604aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 3613b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 362519d70e2SJed Brown ierr = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr); 363084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr); 3643b224e63SBarry Smith 3653b224e63SBarry Smith } else { 3661cee3971SBarry Smith KSP ksp; 3676c924f48SJed Brown char schurprefix[256]; 3683b224e63SBarry Smith 3693b224e63SBarry Smith /* extract the B and C matrices */ 3703b224e63SBarry Smith ilink = jac->head; 3714aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 3724aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 3733b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 3743b224e63SBarry Smith ilink = ilink->next; 3754aa3045dSJed Brown ierr = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr); 3764aa3045dSJed Brown ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 3773b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 378084e4875SJed Brown /* Better would be to use 'mat[0]' (diagonal block of the real matrix) preconditioned by pmat[0] */ 379519d70e2SJed Brown ierr = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr); 3801cee3971SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 381e69d4d44SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr); 3823b224e63SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 3833b224e63SBarry Smith 3843b224e63SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr); 3851cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 386084e4875SJed Brown ierr = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 387084e4875SJed Brown if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) { 388084e4875SJed Brown PC pc; 389cd5f4a64SJed Brown ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr); 390084e4875SJed Brown ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 391084e4875SJed Brown /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */ 392e69d4d44SBarry Smith } 3936c924f48SJed Brown ierr = PetscSNPrintf(schurprefix,sizeof schurprefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",ilink->splitname);CHKERRQ(ierr); 3946c924f48SJed Brown ierr = KSPSetOptionsPrefix(jac->kspschur,schurprefix);CHKERRQ(ierr); 3953b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 3963b224e63SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 3973b224e63SBarry Smith 3983b224e63SBarry Smith ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr); 3993b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr); 4003b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr); 4013b224e63SBarry Smith ilink = jac->head; 4023b224e63SBarry Smith ilink->x = jac->x[0]; ilink->y = jac->y[0]; 4033b224e63SBarry Smith ilink = ilink->next; 4043b224e63SBarry Smith ilink->x = jac->x[1]; ilink->y = jac->y[1]; 4053b224e63SBarry Smith } 4063b224e63SBarry Smith } else { 40797bbdb24SBarry Smith /* set up the individual PCs */ 40897bbdb24SBarry Smith i = 0; 4095a9f2f41SSatish Balay ilink = jac->head; 4105a9f2f41SSatish Balay while (ilink) { 411519d70e2SJed Brown ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr); 4123b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 4135a9f2f41SSatish Balay ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr); 4145a9f2f41SSatish Balay ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 41597bbdb24SBarry Smith i++; 4165a9f2f41SSatish Balay ilink = ilink->next; 4170971522cSBarry Smith } 41897bbdb24SBarry Smith 41997bbdb24SBarry Smith /* create work vectors for each split */ 4201b9fc7fcSBarry Smith if (!jac->x) { 42197bbdb24SBarry Smith ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 4225a9f2f41SSatish Balay ilink = jac->head; 42397bbdb24SBarry Smith for (i=0; i<nsplit; i++) { 424906ed7ccSBarry Smith Vec *vl,*vr; 425906ed7ccSBarry Smith 426906ed7ccSBarry Smith ierr = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr); 427906ed7ccSBarry Smith ilink->x = *vr; 428906ed7ccSBarry Smith ilink->y = *vl; 429906ed7ccSBarry Smith ierr = PetscFree(vr);CHKERRQ(ierr); 430906ed7ccSBarry Smith ierr = PetscFree(vl);CHKERRQ(ierr); 4315a9f2f41SSatish Balay jac->x[i] = ilink->x; 4325a9f2f41SSatish Balay jac->y[i] = ilink->y; 4335a9f2f41SSatish Balay ilink = ilink->next; 43497bbdb24SBarry Smith } 4353b224e63SBarry Smith } 4363b224e63SBarry Smith } 4373b224e63SBarry Smith 4383b224e63SBarry Smith 439d0f46423SBarry Smith if (!jac->head->sctx) { 4403b224e63SBarry Smith Vec xtmp; 4413b224e63SBarry Smith 44279416396SBarry Smith /* compute scatter contexts needed by multiplicative versions and non-default splits */ 4431b9fc7fcSBarry Smith 4445a9f2f41SSatish Balay ilink = jac->head; 4451b9fc7fcSBarry Smith ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 4461b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 447704ba839SBarry Smith ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 4485a9f2f41SSatish Balay ilink = ilink->next; 4491b9fc7fcSBarry Smith } 4501b9fc7fcSBarry Smith ierr = VecDestroy(xtmp);CHKERRQ(ierr); 4511b9fc7fcSBarry Smith } 4520971522cSBarry Smith PetscFunctionReturn(0); 4530971522cSBarry Smith } 4540971522cSBarry Smith 4555a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 456ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 457ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 4585a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 459ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 460ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 46179416396SBarry Smith 4620971522cSBarry Smith #undef __FUNCT__ 4633b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 4643b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 4653b224e63SBarry Smith { 4663b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 4673b224e63SBarry Smith PetscErrorCode ierr; 4683b224e63SBarry Smith KSP ksp; 4693b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 4703b224e63SBarry Smith 4713b224e63SBarry Smith PetscFunctionBegin; 4723b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 4733b224e63SBarry Smith 4743b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4753b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4763b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 4773b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 4783b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 4793b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4803b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4813b224e63SBarry Smith 4823b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 4833b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4843b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4853b224e63SBarry Smith 4863b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 4873b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 4883b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 4893b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4903b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4913b224e63SBarry Smith 4923b224e63SBarry Smith PetscFunctionReturn(0); 4933b224e63SBarry Smith } 4943b224e63SBarry Smith 4953b224e63SBarry Smith #undef __FUNCT__ 4960971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 4970971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 4980971522cSBarry Smith { 4990971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 5000971522cSBarry Smith PetscErrorCode ierr; 5015a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 502af93645bSJed Brown PetscInt cnt; 5030971522cSBarry Smith 5040971522cSBarry Smith PetscFunctionBegin; 50551f519a2SBarry Smith CHKMEMQ; 50651f519a2SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 50751f519a2SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 50851f519a2SBarry Smith 50979416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 5101b9fc7fcSBarry Smith if (jac->defaultsplit) { 5110971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 5125a9f2f41SSatish Balay while (ilink) { 5135a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 5145a9f2f41SSatish Balay ilink = ilink->next; 5150971522cSBarry Smith } 5160971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 5171b9fc7fcSBarry Smith } else { 518efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 5195a9f2f41SSatish Balay while (ilink) { 5205a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 5215a9f2f41SSatish Balay ilink = ilink->next; 5221b9fc7fcSBarry Smith } 5231b9fc7fcSBarry Smith } 52416913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 52579416396SBarry Smith if (!jac->w1) { 52679416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 52779416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 52879416396SBarry Smith } 529efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 5305a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 5313e197d65SBarry Smith cnt = 1; 5325a9f2f41SSatish Balay while (ilink->next) { 5335a9f2f41SSatish Balay ilink = ilink->next; 5343e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 5353e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 5363e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 5373e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5383e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5393e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 5403e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5413e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5423e197d65SBarry Smith } 54351f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 54411755939SBarry Smith cnt -= 2; 54551f519a2SBarry Smith while (ilink->previous) { 54651f519a2SBarry Smith ilink = ilink->previous; 54711755939SBarry Smith /* compute the residual only over the part of the vector needed */ 54811755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 54911755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 55011755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 55111755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 55211755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 55311755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 55411755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 55551f519a2SBarry Smith } 55611755939SBarry Smith } 55765e19b50SBarry Smith } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 55851f519a2SBarry Smith CHKMEMQ; 5590971522cSBarry Smith PetscFunctionReturn(0); 5600971522cSBarry Smith } 5610971522cSBarry Smith 562421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 563ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 564ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 565421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 566ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 567ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 568421e10b8SBarry Smith 569421e10b8SBarry Smith #undef __FUNCT__ 570421e10b8SBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 571421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 572421e10b8SBarry Smith { 573421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 574421e10b8SBarry Smith PetscErrorCode ierr; 575421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 576421e10b8SBarry Smith 577421e10b8SBarry Smith PetscFunctionBegin; 578421e10b8SBarry Smith CHKMEMQ; 579421e10b8SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 580421e10b8SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 581421e10b8SBarry Smith 582421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 583421e10b8SBarry Smith if (jac->defaultsplit) { 584421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 585421e10b8SBarry Smith while (ilink) { 586421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 587421e10b8SBarry Smith ilink = ilink->next; 588421e10b8SBarry Smith } 589421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 590421e10b8SBarry Smith } else { 591421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 592421e10b8SBarry Smith while (ilink) { 593421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 594421e10b8SBarry Smith ilink = ilink->next; 595421e10b8SBarry Smith } 596421e10b8SBarry Smith } 597421e10b8SBarry Smith } else { 598421e10b8SBarry Smith if (!jac->w1) { 599421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 600421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 601421e10b8SBarry Smith } 602421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 603421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 604421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 605421e10b8SBarry Smith while (ilink->next) { 606421e10b8SBarry Smith ilink = ilink->next; 6079989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 608421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 609421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 610421e10b8SBarry Smith } 611421e10b8SBarry Smith while (ilink->previous) { 612421e10b8SBarry Smith ilink = ilink->previous; 6139989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 614421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 615421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 616421e10b8SBarry Smith } 617421e10b8SBarry Smith } else { 618421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 619421e10b8SBarry Smith ilink = ilink->next; 620421e10b8SBarry Smith } 621421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 622421e10b8SBarry Smith while (ilink->previous) { 623421e10b8SBarry Smith ilink = ilink->previous; 6249989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 625421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 626421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 627421e10b8SBarry Smith } 628421e10b8SBarry Smith } 629421e10b8SBarry Smith } 630421e10b8SBarry Smith CHKMEMQ; 631421e10b8SBarry Smith PetscFunctionReturn(0); 632421e10b8SBarry Smith } 633421e10b8SBarry Smith 6340971522cSBarry Smith #undef __FUNCT__ 6350971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 6360971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 6370971522cSBarry Smith { 6380971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6390971522cSBarry Smith PetscErrorCode ierr; 6405a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 6410971522cSBarry Smith 6420971522cSBarry Smith PetscFunctionBegin; 6435a9f2f41SSatish Balay while (ilink) { 6445a9f2f41SSatish Balay ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr); 6455a9f2f41SSatish Balay if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);} 6465a9f2f41SSatish Balay if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);} 6475a9f2f41SSatish Balay if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);} 648704ba839SBarry Smith if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);} 6495a9f2f41SSatish Balay next = ilink->next; 650*7e8c30b6SJed Brown ierr = PetscFree(ilink->splitname);CHKERRQ(ierr); 651704ba839SBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 652704ba839SBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 6535a9f2f41SSatish Balay ilink = next; 6540971522cSBarry Smith } 65505b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 656519d70e2SJed Brown if (jac->mat && jac->mat != jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);} 65797bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 65868dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 65979416396SBarry Smith if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);} 66079416396SBarry Smith if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);} 6613b224e63SBarry Smith if (jac->schur) {ierr = MatDestroy(jac->schur);CHKERRQ(ierr);} 662084e4875SJed Brown if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);} 663d0f46423SBarry Smith if (jac->kspschur) {ierr = KSPDestroy(jac->kspschur);CHKERRQ(ierr);} 6643b224e63SBarry Smith if (jac->B) {ierr = MatDestroy(jac->B);CHKERRQ(ierr);} 6653b224e63SBarry Smith if (jac->C) {ierr = MatDestroy(jac->C);CHKERRQ(ierr);} 6660971522cSBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 6670971522cSBarry Smith PetscFunctionReturn(0); 6680971522cSBarry Smith } 6690971522cSBarry Smith 6700971522cSBarry Smith #undef __FUNCT__ 6710971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 6720971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 6730971522cSBarry Smith { 6741b9fc7fcSBarry Smith PetscErrorCode ierr; 6756c924f48SJed Brown PetscInt bs; 676084e4875SJed Brown PetscTruth flg; 6779dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6783b224e63SBarry Smith PCCompositeType ctype; 6791b9fc7fcSBarry Smith 6800971522cSBarry Smith PetscFunctionBegin; 6811b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 682519d70e2SJed Brown ierr = PetscOptionsTruth("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr); 68351f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 68451f519a2SBarry Smith if (flg) { 68551f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 68651f519a2SBarry Smith } 687704ba839SBarry Smith 6883b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 6893b224e63SBarry Smith if (flg) { 6903b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 6913b224e63SBarry Smith } 692d32f9abdSBarry Smith 693c30613efSMatthew Knepley /* Only setup fields once */ 694c30613efSMatthew Knepley if ((jac->bs > 0) && (jac->nsplits == 0)) { 695d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 696d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 6976c924f48SJed Brown ierr = PCFieldSplitSetRuntimeSplits_Private(pc);CHKERRQ(ierr); 6986c924f48SJed Brown if (jac->splitdefined) {ierr = PetscInfo(pc,"Splits defined using the options database\n");CHKERRQ(ierr);} 699d32f9abdSBarry Smith } 700084e4875SJed 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); 7011b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 7020971522cSBarry Smith PetscFunctionReturn(0); 7030971522cSBarry Smith } 7040971522cSBarry Smith 7050971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 7060971522cSBarry Smith 7070971522cSBarry Smith EXTERN_C_BEGIN 7080971522cSBarry Smith #undef __FUNCT__ 7090971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 7106685144eSJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,const char splitname[],PetscInt n,const PetscInt *fields) 7110971522cSBarry Smith { 71297bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7130971522cSBarry Smith PetscErrorCode ierr; 7145a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 71569a612a9SBarry Smith char prefix[128]; 71651f519a2SBarry Smith PetscInt i; 7170971522cSBarry Smith 7180971522cSBarry Smith PetscFunctionBegin; 7196c924f48SJed Brown if (jac->splitdefined) { 7206c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 7216c924f48SJed Brown PetscFunctionReturn(0); 7226c924f48SJed Brown } 72351f519a2SBarry Smith for (i=0; i<n; i++) { 724e32f2f54SBarry Smith if (fields[i] >= jac->bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs); 725e32f2f54SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 72651f519a2SBarry Smith } 727704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 728db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 729704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 7305a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 7315a9f2f41SSatish Balay ilink->nfields = n; 7325a9f2f41SSatish Balay ilink->next = PETSC_NULL; 7337adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 7341cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 7355a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 73669a612a9SBarry Smith 7376c924f48SJed Brown ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",splitname);CHKERRQ(ierr); 7385a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 7390971522cSBarry Smith 7400971522cSBarry Smith if (!next) { 7415a9f2f41SSatish Balay jac->head = ilink; 74251f519a2SBarry Smith ilink->previous = PETSC_NULL; 7430971522cSBarry Smith } else { 7440971522cSBarry Smith while (next->next) { 7450971522cSBarry Smith next = next->next; 7460971522cSBarry Smith } 7475a9f2f41SSatish Balay next->next = ilink; 74851f519a2SBarry Smith ilink->previous = next; 7490971522cSBarry Smith } 7500971522cSBarry Smith jac->nsplits++; 7510971522cSBarry Smith PetscFunctionReturn(0); 7520971522cSBarry Smith } 7530971522cSBarry Smith EXTERN_C_END 7540971522cSBarry Smith 755e69d4d44SBarry Smith EXTERN_C_BEGIN 756e69d4d44SBarry Smith #undef __FUNCT__ 757e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 758e69d4d44SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 759e69d4d44SBarry Smith { 760e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 761e69d4d44SBarry Smith PetscErrorCode ierr; 762e69d4d44SBarry Smith 763e69d4d44SBarry Smith PetscFunctionBegin; 764519d70e2SJed Brown ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr); 765e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 766e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 767084e4875SJed Brown *n = jac->nsplits; 768e69d4d44SBarry Smith PetscFunctionReturn(0); 769e69d4d44SBarry Smith } 770e69d4d44SBarry Smith EXTERN_C_END 7710971522cSBarry Smith 7720971522cSBarry Smith EXTERN_C_BEGIN 7730971522cSBarry Smith #undef __FUNCT__ 77469a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 775dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 7760971522cSBarry Smith { 7770971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7780971522cSBarry Smith PetscErrorCode ierr; 7790971522cSBarry Smith PetscInt cnt = 0; 7805a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 7810971522cSBarry Smith 7820971522cSBarry Smith PetscFunctionBegin; 78369a612a9SBarry Smith ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 7845a9f2f41SSatish Balay while (ilink) { 7855a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 7865a9f2f41SSatish Balay ilink = ilink->next; 7870971522cSBarry Smith } 788e32f2f54SBarry Smith if (cnt != jac->nsplits) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits); 7890971522cSBarry Smith *n = jac->nsplits; 7900971522cSBarry Smith PetscFunctionReturn(0); 7910971522cSBarry Smith } 7920971522cSBarry Smith EXTERN_C_END 7930971522cSBarry Smith 794704ba839SBarry Smith EXTERN_C_BEGIN 795704ba839SBarry Smith #undef __FUNCT__ 796704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 797db4c96c1SJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,const char splitname[],IS is) 798704ba839SBarry Smith { 799704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 800704ba839SBarry Smith PetscErrorCode ierr; 801704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 802704ba839SBarry Smith char prefix[128]; 803704ba839SBarry Smith 804704ba839SBarry Smith PetscFunctionBegin; 8056c924f48SJed Brown if (jac->splitdefined) { 8066c924f48SJed Brown ierr = PetscInfo1(pc,"Ignoring new split \"%s\" because the splits have already been defined\n",splitname);CHKERRQ(ierr); 8076c924f48SJed Brown PetscFunctionReturn(0); 8086c924f48SJed Brown } 80916913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 810db4c96c1SJed Brown ierr = PetscStrallocpy(splitname,&ilink->splitname);CHKERRQ(ierr); 8111ab39975SBarry Smith ilink->is = is; 8121ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 813704ba839SBarry Smith ilink->next = PETSC_NULL; 814704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 8151cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 816704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 817704ba839SBarry Smith 8186c924f48SJed Brown ierr = PetscSNPrintf(prefix,sizeof prefix,"%sfieldsplit_%s_",((PetscObject)pc)->prefix?((PetscObject)pc)->prefix:"",splitname);CHKERRQ(ierr); 819704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 820704ba839SBarry Smith 821704ba839SBarry Smith if (!next) { 822704ba839SBarry Smith jac->head = ilink; 823704ba839SBarry Smith ilink->previous = PETSC_NULL; 824704ba839SBarry Smith } else { 825704ba839SBarry Smith while (next->next) { 826704ba839SBarry Smith next = next->next; 827704ba839SBarry Smith } 828704ba839SBarry Smith next->next = ilink; 829704ba839SBarry Smith ilink->previous = next; 830704ba839SBarry Smith } 831704ba839SBarry Smith jac->nsplits++; 832704ba839SBarry Smith 833704ba839SBarry Smith PetscFunctionReturn(0); 834704ba839SBarry Smith } 835704ba839SBarry Smith EXTERN_C_END 836704ba839SBarry Smith 8370971522cSBarry Smith #undef __FUNCT__ 8380971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 8390971522cSBarry Smith /*@ 8400971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 8410971522cSBarry Smith 8420971522cSBarry Smith Collective on PC 8430971522cSBarry Smith 8440971522cSBarry Smith Input Parameters: 8450971522cSBarry Smith + pc - the preconditioner context 846db4c96c1SJed Brown . splitname - name of this split 8470971522cSBarry Smith . n - the number of fields in this split 848db4c96c1SJed Brown - fields - the fields in this split 8490971522cSBarry Smith 8500971522cSBarry Smith Level: intermediate 8510971522cSBarry Smith 852d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 853d32f9abdSBarry Smith 854d32f9abdSBarry Smith The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block 855d32f9abdSBarry 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 856d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 857d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 858d32f9abdSBarry Smith 859db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 860db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 861db4c96c1SJed Brown 862d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 8630971522cSBarry Smith 8640971522cSBarry Smith @*/ 8656685144eSJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,const char splitname[],PetscInt n,const PetscInt *fields) 8660971522cSBarry Smith { 8676685144eSJed Brown PetscErrorCode ierr,(*f)(PC,const char[],PetscInt,const PetscInt *); 8680971522cSBarry Smith 8690971522cSBarry Smith PetscFunctionBegin; 8700700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 871db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 872db4c96c1SJed Brown if (n < 1) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Provided number of fields %D in split \"%s\" not positive",n,splitname); 873db4c96c1SJed Brown PetscValidIntPointer(fields,3); 8740971522cSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr); 8750971522cSBarry Smith if (f) { 876db4c96c1SJed Brown ierr = (*f)(pc,splitname,n,fields);CHKERRQ(ierr); 8770971522cSBarry Smith } 8780971522cSBarry Smith PetscFunctionReturn(0); 8790971522cSBarry Smith } 8800971522cSBarry Smith 8810971522cSBarry Smith #undef __FUNCT__ 882704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 883704ba839SBarry Smith /*@ 884704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 885704ba839SBarry Smith 886704ba839SBarry Smith Collective on PC 887704ba839SBarry Smith 888704ba839SBarry Smith Input Parameters: 889704ba839SBarry Smith + pc - the preconditioner context 890db4c96c1SJed Brown . splitname - name of this split 891db4c96c1SJed Brown - is - the index set that defines the vector elements in this field 892704ba839SBarry Smith 893d32f9abdSBarry Smith 894a6ffb8dbSJed Brown Notes: 895a6ffb8dbSJed Brown Use PCFieldSplitSetFields(), for fields defined by strided types. 896a6ffb8dbSJed Brown 897db4c96c1SJed Brown This function is called once per split (it creates a new split each time). Solve options 898db4c96c1SJed Brown for this split will be available under the prefix -fieldsplit_SPLITNAME_. 899d32f9abdSBarry Smith 900704ba839SBarry Smith Level: intermediate 901704ba839SBarry Smith 902704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 903704ba839SBarry Smith 904704ba839SBarry Smith @*/ 905db4c96c1SJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,const char splitname[],IS is) 906704ba839SBarry Smith { 907db4c96c1SJed Brown PetscErrorCode ierr,(*f)(PC,const char[],IS); 908704ba839SBarry Smith 909704ba839SBarry Smith PetscFunctionBegin; 9100700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 911db4c96c1SJed Brown PetscValidCharPointer(splitname,2); 912db4c96c1SJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,3); 913704ba839SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr); 914704ba839SBarry Smith if (f) { 915db4c96c1SJed Brown ierr = (*f)(pc,splitname,is);CHKERRQ(ierr); 916704ba839SBarry Smith } 917704ba839SBarry Smith PetscFunctionReturn(0); 918704ba839SBarry Smith } 919704ba839SBarry Smith 920704ba839SBarry Smith #undef __FUNCT__ 92151f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 92251f519a2SBarry Smith /*@ 92351f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 92451f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 92551f519a2SBarry Smith 92651f519a2SBarry Smith Collective on PC 92751f519a2SBarry Smith 92851f519a2SBarry Smith Input Parameters: 92951f519a2SBarry Smith + pc - the preconditioner context 93051f519a2SBarry Smith - bs - the block size 93151f519a2SBarry Smith 93251f519a2SBarry Smith Level: intermediate 93351f519a2SBarry Smith 93451f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 93551f519a2SBarry Smith 93651f519a2SBarry Smith @*/ 93751f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 93851f519a2SBarry Smith { 93951f519a2SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt); 94051f519a2SBarry Smith 94151f519a2SBarry Smith PetscFunctionBegin; 9420700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 94351f519a2SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr); 94451f519a2SBarry Smith if (f) { 94551f519a2SBarry Smith ierr = (*f)(pc,bs);CHKERRQ(ierr); 94651f519a2SBarry Smith } 94751f519a2SBarry Smith PetscFunctionReturn(0); 94851f519a2SBarry Smith } 94951f519a2SBarry Smith 95051f519a2SBarry Smith #undef __FUNCT__ 95169a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 9520971522cSBarry Smith /*@C 95369a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 9540971522cSBarry Smith 95569a612a9SBarry Smith Collective on KSP 9560971522cSBarry Smith 9570971522cSBarry Smith Input Parameter: 9580971522cSBarry Smith . pc - the preconditioner context 9590971522cSBarry Smith 9600971522cSBarry Smith Output Parameters: 9610971522cSBarry Smith + n - the number of split 96269a612a9SBarry Smith - pc - the array of KSP contexts 9630971522cSBarry Smith 9640971522cSBarry Smith Note: 965d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 966d32f9abdSBarry Smith (not the KSP just the array that contains them). 9670971522cSBarry Smith 96869a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 9690971522cSBarry Smith 9700971522cSBarry Smith Level: advanced 9710971522cSBarry Smith 9720971522cSBarry Smith .seealso: PCFIELDSPLIT 9730971522cSBarry Smith @*/ 974dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 9750971522cSBarry Smith { 97669a612a9SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **); 9770971522cSBarry Smith 9780971522cSBarry Smith PetscFunctionBegin; 9790700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9800971522cSBarry Smith PetscValidIntPointer(n,2); 98169a612a9SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 9820971522cSBarry Smith if (f) { 98369a612a9SBarry Smith ierr = (*f)(pc,n,subksp);CHKERRQ(ierr); 984e7e72b3dSBarry Smith } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC"); 9850971522cSBarry Smith PetscFunctionReturn(0); 9860971522cSBarry Smith } 9870971522cSBarry Smith 988e69d4d44SBarry Smith #undef __FUNCT__ 989e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 990e69d4d44SBarry Smith /*@ 991e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 992e69d4d44SBarry Smith D matrix. Otherwise no preconditioner is used. 993e69d4d44SBarry Smith 994e69d4d44SBarry Smith Collective on PC 995e69d4d44SBarry Smith 996e69d4d44SBarry Smith Input Parameters: 997e69d4d44SBarry Smith + pc - the preconditioner context 998084e4875SJed Brown . ptype - which matrix to use for preconditioning the Schur complement 999084e4875SJed Brown - userpre - matrix to use for preconditioning, or PETSC_NULL 1000084e4875SJed Brown 1001084e4875SJed Brown Notes: 1002084e4875SJed Brown The default is to use the block on the diagonal of the preconditioning matrix. This is D, in the (1,1) position. 1003084e4875SJed Brown There are currently no preconditioners that work directly with the Schur complement so setting 1004084e4875SJed Brown PC_FIELDSPLIT_SCHUR_PRE_SELF is observationally equivalent to -fieldsplit_1_pc_type none. 1005e69d4d44SBarry Smith 1006e69d4d44SBarry Smith Options Database: 1007084e4875SJed Brown . -pc_fieldsplit_schur_precondition <self,user,diag> default is diag 1008e69d4d44SBarry Smith 1009e69d4d44SBarry Smith Level: intermediate 1010e69d4d44SBarry Smith 1011084e4875SJed Brown .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType 1012e69d4d44SBarry Smith 1013e69d4d44SBarry Smith @*/ 1014084e4875SJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1015e69d4d44SBarry Smith { 1016084e4875SJed Brown PetscErrorCode ierr,(*f)(PC,PCFieldSplitSchurPreType,Mat); 1017e69d4d44SBarry Smith 1018e69d4d44SBarry Smith PetscFunctionBegin; 10190700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1020e69d4d44SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",(void (**)(void))&f);CHKERRQ(ierr); 1021e69d4d44SBarry Smith if (f) { 1022084e4875SJed Brown ierr = (*f)(pc,ptype,pre);CHKERRQ(ierr); 1023e69d4d44SBarry Smith } 1024e69d4d44SBarry Smith PetscFunctionReturn(0); 1025e69d4d44SBarry Smith } 1026e69d4d44SBarry Smith 1027e69d4d44SBarry Smith EXTERN_C_BEGIN 1028e69d4d44SBarry Smith #undef __FUNCT__ 1029e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 1030084e4875SJed Brown PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre) 1031e69d4d44SBarry Smith { 1032e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1033084e4875SJed Brown PetscErrorCode ierr; 1034e69d4d44SBarry Smith 1035e69d4d44SBarry Smith PetscFunctionBegin; 1036084e4875SJed Brown jac->schurpre = ptype; 1037084e4875SJed Brown if (pre) { 1038084e4875SJed Brown if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);} 1039084e4875SJed Brown jac->schur_user = pre; 1040084e4875SJed Brown ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr); 1041084e4875SJed Brown } 1042e69d4d44SBarry Smith PetscFunctionReturn(0); 1043e69d4d44SBarry Smith } 1044e69d4d44SBarry Smith EXTERN_C_END 1045e69d4d44SBarry Smith 104630ad9308SMatthew Knepley #undef __FUNCT__ 104730ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 104830ad9308SMatthew Knepley /*@C 104930ad9308SMatthew Knepley PCFieldSplitGetSchurBlocks - Gets the all matrix blocks for the Schur complement 105030ad9308SMatthew Knepley 105130ad9308SMatthew Knepley Collective on KSP 105230ad9308SMatthew Knepley 105330ad9308SMatthew Knepley Input Parameter: 105430ad9308SMatthew Knepley . pc - the preconditioner context 105530ad9308SMatthew Knepley 105630ad9308SMatthew Knepley Output Parameters: 105730ad9308SMatthew Knepley + A - the (0,0) block 105830ad9308SMatthew Knepley . B - the (0,1) block 105930ad9308SMatthew Knepley . C - the (1,0) block 106030ad9308SMatthew Knepley - D - the (1,1) block 106130ad9308SMatthew Knepley 106230ad9308SMatthew Knepley Level: advanced 106330ad9308SMatthew Knepley 106430ad9308SMatthew Knepley .seealso: PCFIELDSPLIT 106530ad9308SMatthew Knepley @*/ 106630ad9308SMatthew Knepley PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSchurBlocks(PC pc,Mat *A,Mat *B,Mat *C, Mat *D) 106730ad9308SMatthew Knepley { 106830ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 106930ad9308SMatthew Knepley 107030ad9308SMatthew Knepley PetscFunctionBegin; 10710700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1072c1235816SBarry Smith if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach."); 107330ad9308SMatthew Knepley if (A) *A = jac->pmat[0]; 107430ad9308SMatthew Knepley if (B) *B = jac->B; 107530ad9308SMatthew Knepley if (C) *C = jac->C; 107630ad9308SMatthew Knepley if (D) *D = jac->pmat[1]; 107730ad9308SMatthew Knepley PetscFunctionReturn(0); 107830ad9308SMatthew Knepley } 107930ad9308SMatthew Knepley 108079416396SBarry Smith EXTERN_C_BEGIN 108179416396SBarry Smith #undef __FUNCT__ 108279416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 1083dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 108479416396SBarry Smith { 108579416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1086e69d4d44SBarry Smith PetscErrorCode ierr; 108779416396SBarry Smith 108879416396SBarry Smith PetscFunctionBegin; 108979416396SBarry Smith jac->type = type; 10903b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 10913b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 10923b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 1093e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1094e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1095e69d4d44SBarry Smith 10963b224e63SBarry Smith } else { 10973b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 10983b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 1099e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 11009e7d6b0aSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 11013b224e63SBarry Smith } 110279416396SBarry Smith PetscFunctionReturn(0); 110379416396SBarry Smith } 110479416396SBarry Smith EXTERN_C_END 110579416396SBarry Smith 110651f519a2SBarry Smith EXTERN_C_BEGIN 110751f519a2SBarry Smith #undef __FUNCT__ 110851f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 110951f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 111051f519a2SBarry Smith { 111151f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 111251f519a2SBarry Smith 111351f519a2SBarry Smith PetscFunctionBegin; 111465e19b50SBarry Smith if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 111565e19b50SBarry Smith if (jac->bs > 0 && jac->bs != bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot change fieldsplit blocksize from %D to %D after it has been set",jac->bs,bs); 111651f519a2SBarry Smith jac->bs = bs; 111751f519a2SBarry Smith PetscFunctionReturn(0); 111851f519a2SBarry Smith } 111951f519a2SBarry Smith EXTERN_C_END 112051f519a2SBarry Smith 112179416396SBarry Smith #undef __FUNCT__ 112279416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1123bc08b0f1SBarry Smith /*@ 112479416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 112579416396SBarry Smith 112679416396SBarry Smith Collective on PC 112779416396SBarry Smith 112879416396SBarry Smith Input Parameter: 112979416396SBarry Smith . pc - the preconditioner context 113081540f2fSBarry Smith . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 113179416396SBarry Smith 113279416396SBarry Smith Options Database Key: 1133a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 113479416396SBarry Smith 113579416396SBarry Smith Level: Developer 113679416396SBarry Smith 113779416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 113879416396SBarry Smith 113979416396SBarry Smith .seealso: PCCompositeSetType() 114079416396SBarry Smith 114179416396SBarry Smith @*/ 1142dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type) 114379416396SBarry Smith { 114479416396SBarry Smith PetscErrorCode ierr,(*f)(PC,PCCompositeType); 114579416396SBarry Smith 114679416396SBarry Smith PetscFunctionBegin; 11470700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 114879416396SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 114979416396SBarry Smith if (f) { 115079416396SBarry Smith ierr = (*f)(pc,type);CHKERRQ(ierr); 115179416396SBarry Smith } 115279416396SBarry Smith PetscFunctionReturn(0); 115379416396SBarry Smith } 115479416396SBarry Smith 11550971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 11560971522cSBarry Smith /*MC 1157a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 11580971522cSBarry Smith fields or groups of fields 11590971522cSBarry Smith 11600971522cSBarry Smith 1161edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1162edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 11630971522cSBarry Smith 1164a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 116569a612a9SBarry Smith and set the options directly on the resulting KSP object 11660971522cSBarry Smith 11670971522cSBarry Smith Level: intermediate 11680971522cSBarry Smith 116979416396SBarry Smith Options Database Keys: 117081540f2fSBarry Smith + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 117181540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 117281540f2fSBarry Smith been supplied explicitly by -pc_fieldsplit_%d_fields 117381540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 117481540f2fSBarry Smith . -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative> 1175e69d4d44SBarry Smith . -pc_fieldsplit_schur_precondition <true,false> default is true 117679416396SBarry Smith 1177edf189efSBarry Smith - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 11783b224e63SBarry Smith for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 11793b224e63SBarry Smith 11803b224e63SBarry Smith 1181d32f9abdSBarry Smith Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1182d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1183d32f9abdSBarry Smith 1184d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1185d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1186d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1187d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1188d32f9abdSBarry Smith 1189d32f9abdSBarry Smith Currently for the multiplicative version, the updated residual needed for the next field 1190d32f9abdSBarry Smith solve is computed via a matrix vector product over the entire array. An optimization would be 1191d32f9abdSBarry Smith to update the residual only for the part of the right hand side associated with the next field 1192d32f9abdSBarry Smith solve. (This would involve more MatGetSubMatrix() calls or some other mechanism to compute the 1193d32f9abdSBarry Smith part of the matrix needed to just update part of the residual). 1194d32f9abdSBarry Smith 1195e69d4d44SBarry Smith For the Schur complement preconditioner if J = ( A B ) 1196e69d4d44SBarry Smith ( C D ) 1197e69d4d44SBarry Smith the preconditioner is 1198e69d4d44SBarry Smith (I -B inv(A)) ( inv(A) 0 ) (I 0 ) 1199e69d4d44SBarry Smith (0 I ) ( 0 inv(S) ) (-C inv(A) I ) 1200edf189efSBarry Smith where the action of inv(A) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1201e69d4d44SBarry Smith inv(S) is computed using the KSP solver with prefix -schur_. For PCFieldSplitGetKSP() when field number is 1202edf189efSBarry Smith 0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1203e69d4d44SBarry Smith D is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1204e69d4d44SBarry Smith option. 1205e69d4d44SBarry Smith 1206edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1207edf189efSBarry Smith is used automatically for a second block. 1208edf189efSBarry Smith 1209a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 12100971522cSBarry Smith 1211a541d17aSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners 1212e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 12130971522cSBarry Smith M*/ 12140971522cSBarry Smith 12150971522cSBarry Smith EXTERN_C_BEGIN 12160971522cSBarry Smith #undef __FUNCT__ 12170971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 1218dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc) 12190971522cSBarry Smith { 12200971522cSBarry Smith PetscErrorCode ierr; 12210971522cSBarry Smith PC_FieldSplit *jac; 12220971522cSBarry Smith 12230971522cSBarry Smith PetscFunctionBegin; 122438f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 12250971522cSBarry Smith jac->bs = -1; 12260971522cSBarry Smith jac->nsplits = 0; 12273e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1228e6cab6aaSJed Brown jac->schurpre = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */ 122951f519a2SBarry Smith 12300971522cSBarry Smith pc->data = (void*)jac; 12310971522cSBarry Smith 12320971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1233421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 12340971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 12350971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 12360971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 12370971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 12380971522cSBarry Smith pc->ops->applyrichardson = 0; 12390971522cSBarry Smith 124069a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 124169a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 12420971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 12430971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1244704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1245704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 124679416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 124779416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 124851f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 124951f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 12500971522cSBarry Smith PetscFunctionReturn(0); 12510971522cSBarry Smith } 12520971522cSBarry Smith EXTERN_C_END 12530971522cSBarry Smith 12540971522cSBarry Smith 1255a541d17aSBarry Smith /*MC 1256a541d17aSBarry Smith Block_Preconditioners - PETSc provides support for a variety of "block preconditioners", this provides an 1257a541d17aSBarry Smith overview of these methods. 1258a541d17aSBarry Smith 1259a541d17aSBarry Smith Consider the solution to ( A B ) (x_1) = (b_1) 1260a541d17aSBarry Smith ( C D ) (x_2) (b_2) 1261a541d17aSBarry Smith 1262a541d17aSBarry Smith Important special cases, the Stokes equation: C = B' and D = 0 (A B) (x_1) = (b_1) 1263a541d17aSBarry Smith B' 0) (x_2) (b_2) 1264a541d17aSBarry Smith 1265a541d17aSBarry Smith One of the goals of the PCFieldSplit preconditioner in PETSc is to provide a variety of preconditioners 1266a541d17aSBarry Smith for this block system. 1267a541d17aSBarry Smith 1268a541d17aSBarry Smith Consider an additional matrix (Ap Bp) 1269a541d17aSBarry Smith (Cp Dp) where some or all of the entries may be the same as 1270a541d17aSBarry Smith in the original matrix (for example Ap == A). 1271a541d17aSBarry Smith 1272a541d17aSBarry Smith In the following, A^ denotes the approximate application of the inverse of A, possibly using Ap in the 1273a541d17aSBarry Smith approximation. In PETSc this simply means one has called KSPSetOperators(ksp,A,Ap,...) or KSPSetOperators(ksp,Ap,Ap,...) 1274a541d17aSBarry Smith 1275a541d17aSBarry Smith Block Jacobi: x_1 = A^ b_1 1276a541d17aSBarry Smith x_2 = D^ b_2 1277a541d17aSBarry Smith 1278a541d17aSBarry Smith Lower block Gauss-Seidel: x_1 = A^ b_1 1279a541d17aSBarry Smith x_2 = D^ (b_2 - C x_1) variant x_2 = D^ (b_2 - Cp x_1) 1280a541d17aSBarry Smith 1281a541d17aSBarry 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) 1282a541d17aSBarry Smith Interestingly this form is not actually a symmetric matrix, the symmetric version is 1283a541d17aSBarry Smith x_1 = A^(b_1 - B x_2) variant x_1 = A^(b_1 - Bp x_2) 1284a541d17aSBarry Smith 12850bc0a719SBarry Smith Level: intermediate 12860bc0a719SBarry Smith 1287a541d17aSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCFIELDSPLIT 1288a541d17aSBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 1289a541d17aSBarry Smith M*/ 1290